/[MITgcm]/MITgcm/pkg/exf/exf_mapfields.F
ViewVC logotype

Contents of /MITgcm/pkg/exf/exf_mapfields.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.4 - (show annotations) (download)
Sat Dec 28 10:11:11 2002 UTC (21 years, 5 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint48e_post, checkpoint48b_post, checkpoint48c_pre, checkpoint48d_pre, checkpoint47i_post, checkpoint48d_post, checkpoint47g_post, checkpoint48a_post, checkpoint47j_post, checkpoint48c_post, checkpoint47f_post, checkpoint48, checkpoint47h_post
Changes since 1.3: +5 -5 lines
checkpoint47f_post
Merging from release1_p10:
o modifications for using pkg/exf with pkg/seaice
  - pkg/seaice CPP options SEAICE_EXTERNAL_FORCING
    and SEAICE_EXTERNAL_FLUXES
  - pkg/exf CPP options EXF_READ_EVAP and
    EXF_NO_BULK_COMPUTATIONS
  - usage examples are Experiments 8 and 9 in
    verification/lab_sea/README
  - verification/lab_sea default experiment now uses
    pkg/gmredi, pkg/kpp, pkg/seaice, and pkg/exf

1 c $Header: /u/gcmpack/MITgcm/pkg/exf/exf_mapfields.F,v 1.3 2002/12/17 19:47:41 cheisey Exp $
2
3 #include "EXF_CPPOPTIONS.h"
4
5
6 subroutine exf_mapfields( mythid )
7
8 c ==================================================================
9 c SUBROUTINE exf_mapfields
10 c ==================================================================
11 c
12 c o Map the external forcing fields on the ocean model arrays. This
13 c routine is included to separate the ocean state estimation tool
14 c as much as possible from the ocean model. Unit conversion factors
15 c are to be set by the user.
16 c
17 c The units have to be such that the individual forcing record has
18 c units equal to [quantity/s]. See the header file *FFIELDS.h* of
19 c the MITgcmuv.
20 c
21 c Required units such that no scaling has to be applied:
22 c
23 c heat flux: input file W/m^2
24 c salt flux: input file m/s
25 c zonal wind stress: input file N/m^2
26 c merid. wind stress: input file N/m^2
27 c
28 c To allow for such unit conversions this routine contains scaling
29 c factors scal_quantity.
30 c
31 c started: Christian Eckert eckert@mit.edu 09-Aug-1999
32 c
33 c changed: Christian Eckert eckert@mit.edu 11-Jan-2000
34 c
35 c - Restructured the code in order to create a package
36 c for the MITgcmUV.
37 c
38 c Christian Eckert eckert@mit.edu 12-Feb-2000
39 c
40 c - Changed Routine names (package prefix: exf_)
41 c
42 c Patrick Heimbach, heimbach@mit.edu 06-May-2000
43 c
44 c - added and changed CPP flag structure for
45 c ALLOW_BULKFORMULAE, ALLOW_ATM_TEMP
46 c
47 c Patrick Heimbach, heimbach@mit.edu 23-May-2000
48 c
49 c - sign change of ustress/vstress incorporated into
50 c scaling factors scal_ust, scal_vst
51 c
52 c Dimitris Menemenlis, menemenlis@jpl.nasa.gov 20-Dec-2002
53 c
54 c - removed: empmr(i,j,bi,bj) = scal_prc*precip(i,j,bi,bj)
55 c
56 c ==================================================================
57 c SUBROUTINE exf_mapfields
58 c ==================================================================
59
60 implicit none
61
62 c == global variables ==
63
64 #include "EEPARAMS.h"
65 #include "SIZE.h"
66 #include "FFIELDS.h"
67 #include "exf_param.h"
68 #include "exf_constants.h"
69 #include "exf_fields.h"
70 #include "exf_clim_fields.h"
71 #ifdef ALLOW_AUTODIFF_TAMC
72 # include "tamc.h"
73 # include "tamc_keys.h"
74 #endif
75 c == routine arguments ==
76
77 c mythid - thread number for this instance of the routine.
78
79 integer mythid
80
81 c == local variables ==
82
83 integer bi,bj
84 integer i,j
85 integer jtlo
86 integer jthi
87 integer itlo
88 integer ithi
89 integer jmin
90 integer jmax
91 integer imin
92 integer imax
93
94 c == end of interface ==
95
96 jtlo = mybylo(mythid)
97 jthi = mybyhi(mythid)
98 itlo = mybxlo(mythid)
99 ithi = mybxhi(mythid)
100 jmin = 1-oly
101 jmax = sny+oly
102 imin = 1-olx
103 imax = snx+olx
104
105 do bj = jtlo,jthi
106 do bi = itlo,ithi
107
108 #ifdef ALLOW_AUTODIFF_TAMC
109 act1 = bi - myBxLo(myThid)
110 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
111 act2 = bj - myByLo(myThid)
112 max2 = myByHi(myThid) - myByLo(myThid) + 1
113 act3 = myThid - 1
114 max3 = nTx*nTy
115 act4 = ikey_dynamics - 1
116 ikey = (act1 + 1) + act2*max1
117 & + act3*max1*max2
118 & + act4*max1*max2*max3
119 #endif /* ALLOW_AUTODIFF_TAMC */
120
121 do j = jmin,jmax
122 do i = imin,imax
123 c Heat flux.
124 qnet(i,j,bi,bj) = scal_hfl*hflux(i,j,bi,bj)
125 enddo
126 enddo
127
128
129 do j = jmin,jmax
130 do i = imin,imax
131 c Salt flux.
132 empmr(i,j,bi,bj) = scal_sfl*sflux(i,j,bi,bj)
133 enddo
134 enddo
135
136 #ifdef ALLOW_AUTODIFF_TAMC
137 CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
138 #endif
139 do j = jmin,jmax
140 do i = imin,imax
141 c Zonal wind stress.
142 if (ustress(i,j,bi,bj).gt.2.0D0) then
143 ustress(i,j,bi,bj)=2.0D0
144 endif
145 enddo
146 enddo
147 #ifdef ALLOW_AUTODIFF_TAMC
148 CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
149 #endif
150 do j = jmin,jmax
151 do i = imin,imax
152 if (ustress(i,j,bi,bj).lt.-2.0D0) then
153 ustress(i,j,bi,bj)=-2.0D0
154 endif
155 enddo
156 enddo
157 do j = jmin,jmax
158 do i = imin,imax
159 fu(i,j,bi,bj) = scal_ust*ustress(i,j,bi,bj)
160 enddo
161 enddo
162
163 #ifdef ALLOW_AUTODIFF_TAMC
164 CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
165 #endif
166 do j = jmin,jmax
167 do i = imin,imax
168 c Meridional wind stress.
169 if (vstress(i,j,bi,bj).gt.2.0D0) then
170 vstress(i,j,bi,bj)=2.0D0
171 endif
172 enddo
173 enddo
174 #ifdef ALLOW_AUTODIFF_TAMC
175 CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
176 #endif
177 do j = jmin,jmax
178 do i = imin,imax
179 if (vstress(i,j,bi,bj).lt.-2.0D0) then
180 vstress(i,j,bi,bj)=-2.0D0
181 endif
182 enddo
183 enddo
184 do j = jmin,jmax
185 do i = imin,imax
186 fv(i,j,bi,bj) = scal_vst*vstress(i,j,bi,bj)
187 enddo
188 enddo
189
190 #ifdef ALLOW_KPP || \
191 (defined (ALLOW_BULKFORMULAE) && defined (ALLOW_ATM_TEMP)))
192 c Short wave radiative flux.
193 do j = jmin,jmax
194 do i = imin,imax
195 qsw(i,j,bi,bj) = scal_swf*swflux(i,j,bi,bj)
196 enddo
197 enddo
198 #endif
199
200 #ifdef ALLOW_CLIMSST_RELAXATION
201 do j = jmin,jmax
202 do i = imin,imax
203 sst(i,j,bi,bj) = scal_sst*climsst(i,j,bi,bj)
204 enddo
205 enddo
206 #endif
207
208 #ifdef ALLOW_CLIMSSS_RELAXATION
209 do j = jmin,jmax
210 do i = imin,imax
211 sss(i,j,bi,bj) = scal_sss*climsss(i,j,bi,bj)
212 enddo
213 enddo
214 #endif
215
216 #ifdef ATMOSPHERIC_LOADING
217 do j = jmin,jmax
218 do i = imin,imax
219 pload(i,j,bi,bj) = scal_apressure*apressure(i,j,bi,bj)
220 enddo
221 enddo
222 #endif
223
224
225 enddo
226 enddo
227
228 c Update the tile edges.
229
230 _EXCH_XY_R4( qnet, mythid )
231 _EXCH_XY_R4( empmr, mythid )
232 c _EXCH_XY_R4( fu, mythid )
233 c _EXCH_XY_R4( fv, mythid )
234 CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid)
235 #ifdef ALLOW_KPP
236 _EXCH_XY_R4( qsw, mythid )
237 #endif
238 #ifdef ALLOW_CLIMSST_RELAXATION
239 _EXCH_XY_R4( sst, mythid )
240 #endif
241 #ifdef ALLOW_CLIMSSS_RELAXATION
242 _EXCH_XY_R4( sss, mythid )
243 #endif
244 #ifdef ATMOSPHERIC_LOADING
245 _EXCH_XY_R4( pload, mythid )
246 #endif
247
248 end

  ViewVC Help
Powered by ViewVC 1.1.22