1 |
C $Header: /u/gcmpack/MITgcm/pkg/exf/exf_mapfields.F,v 1.24 2011/12/22 19:03:41 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "EXF_OPTIONS.h" |
5 |
|
6 |
subroutine exf_mapfields( mytime, myiter, mythid ) |
7 |
|
8 |
c ================================================================== |
9 |
c SUBROUTINE exf_mapfields |
10 |
c ================================================================== |
11 |
c |
12 |
c o Map external forcing fields (ustress, vstress, hflux, sflux, |
13 |
c swflux, apressure, climsss, climsst, etc.) onto ocean model |
14 |
c arrays (fu, fv, Qnet, EmPmR, Qsw, pload, sss, sst, etc.). |
15 |
c This routine is included to separate the ocean state estimation |
16 |
c tool as much as possible from the ocean model. Unit and sign |
17 |
c conventions can be customized using variables exf_outscal_*, |
18 |
c which are set in exf_readparms.F. See the header files |
19 |
c EXF_FIELDS.h and FFIELDS.h for definitions of the various input |
20 |
c and output fields and for default unit and sign convetions. |
21 |
c |
22 |
c started: Christian Eckert eckert@mit.edu 09-Aug-1999 |
23 |
c |
24 |
c changed: Christian Eckert eckert@mit.edu 11-Jan-2000 |
25 |
c - Restructured the code in order to create a package |
26 |
c for the MITgcmUV. |
27 |
c |
28 |
c Christian Eckert eckert@mit.edu 12-Feb-2000 |
29 |
c - Changed Routine names (package prefix: exf_) |
30 |
c |
31 |
c Patrick Heimbach, heimbach@mit.edu 06-May-2000 |
32 |
c - added and changed CPP flag structure for |
33 |
c ALLOW_BULKFORMULAE, ALLOW_ATM_TEMP |
34 |
c |
35 |
c Patrick Heimbach, heimbach@mit.edu 23-May-2000 |
36 |
c - sign change of ustress/vstress incorporated into |
37 |
c scaling factors exf_outscal_ust, exf_outscal_vst |
38 |
c |
39 |
c mods for pkg/seaice: menemenlis@jpl.nasa.gov 20-Dec-2002 |
40 |
c |
41 |
c ================================================================== |
42 |
c SUBROUTINE exf_mapfields |
43 |
c ================================================================== |
44 |
|
45 |
implicit none |
46 |
|
47 |
c == global variables == |
48 |
|
49 |
#include "EEPARAMS.h" |
50 |
#include "SIZE.h" |
51 |
#include "PARAMS.h" |
52 |
#include "FFIELDS.h" |
53 |
#include "GRID.h" |
54 |
|
55 |
#include "EXF_PARAM.h" |
56 |
#include "EXF_CONSTANTS.h" |
57 |
#include "EXF_FIELDS.h" |
58 |
#ifdef ALLOW_AUTODIFF_TAMC |
59 |
# include "tamc.h" |
60 |
# include "tamc_keys.h" |
61 |
#endif |
62 |
c == routine arguments == |
63 |
|
64 |
c mythid - thread number for this instance of the routine. |
65 |
|
66 |
_RL mytime |
67 |
integer myiter |
68 |
integer mythid |
69 |
|
70 |
c == local variables == |
71 |
|
72 |
INTEGER bi,bj |
73 |
INTEGER i,j,ks |
74 |
INTEGER imin, imax |
75 |
INTEGER jmin, jmax |
76 |
PARAMETER ( imin = 1-OLx , imax = sNx+OLx ) |
77 |
PARAMETER ( jmin = 1-OLy , jmax = sNy+OLy ) |
78 |
|
79 |
c == end of interface == |
80 |
|
81 |
C-- set surface level index: |
82 |
ks = 1 |
83 |
|
84 |
DO bj = myByLo(myThid),myByHi(myThid) |
85 |
DO bi = myBxLo(myThid),myBxHi(myThid) |
86 |
|
87 |
#ifdef ALLOW_AUTODIFF_TAMC |
88 |
act1 = bi - myBxLo(myThid) |
89 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
90 |
act2 = bj - myByLo(myThid) |
91 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
92 |
act3 = myThid - 1 |
93 |
max3 = nTx*nTy |
94 |
act4 = ikey_dynamics - 1 |
95 |
ikey = (act1 + 1) + act2*max1 |
96 |
& + act3*max1*max2 |
97 |
& + act4*max1*max2*max3 |
98 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
99 |
|
100 |
c Heat flux. |
101 |
do j = jmin,jmax |
102 |
do i = imin,imax |
103 |
qnet(i,j,bi,bj) = exf_outscal_hflux*hflux(i,j,bi,bj) |
104 |
enddo |
105 |
enddo |
106 |
if ( hfluxfile .EQ. ' ' ) then |
107 |
do j = jmin,jmax |
108 |
do i = imin,imax |
109 |
qnet(i,j,bi,bj) = qnet(i,j,bi,bj) - |
110 |
& exf_outscal_hflux * ( hflux_exfremo_intercept + |
111 |
& hflux_exfremo_slope*(mytime-starttime) ) |
112 |
enddo |
113 |
enddo |
114 |
endif |
115 |
|
116 |
c Salt flux. |
117 |
do j = jmin,jmax |
118 |
do i = imin,imax |
119 |
EmPmR(i,j,bi,bj)= exf_outscal_sflux*sflux(i,j,bi,bj) |
120 |
& *rhoConstFresh |
121 |
enddo |
122 |
enddo |
123 |
if ( sfluxfile .EQ. ' ' ) then |
124 |
do j = jmin,jmax |
125 |
do i = imin,imax |
126 |
EmPmR(i,j,bi,bj) = EmPmR(i,j,bi,bj) - rhoConstFresh* |
127 |
& exf_outscal_sflux * ( sflux_exfremo_intercept + |
128 |
& sflux_exfremo_slope*(mytime-starttime) ) |
129 |
enddo |
130 |
enddo |
131 |
endif |
132 |
|
133 |
#ifdef ALLOW_AUTODIFF_TAMC |
134 |
CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
135 |
#endif |
136 |
do j = jmin,jmax |
137 |
do i = imin,imax |
138 |
c Zonal wind stress. |
139 |
if (ustress(i,j,bi,bj).gt.windstressmax) then |
140 |
ustress(i,j,bi,bj)=windstressmax |
141 |
endif |
142 |
enddo |
143 |
enddo |
144 |
#ifdef ALLOW_AUTODIFF_TAMC |
145 |
CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
146 |
#endif |
147 |
do j = jmin,jmax |
148 |
do i = imin,imax |
149 |
if (ustress(i,j,bi,bj).lt.-windstressmax) then |
150 |
ustress(i,j,bi,bj)=-windstressmax |
151 |
endif |
152 |
enddo |
153 |
enddo |
154 |
IF ( stressIsOnCgrid ) THEN |
155 |
do j = jmin,jmax |
156 |
do i = imin+1,imax |
157 |
fu(i,j,bi,bj) = exf_outscal_ustress*ustress(i,j,bi,bj) |
158 |
enddo |
159 |
enddo |
160 |
ELSE |
161 |
do j = jmin,jmax |
162 |
do i = imin+1,imax |
163 |
c Shift wind stresses calculated at Grid-center to W/S points |
164 |
fu(i,j,bi,bj) = exf_outscal_ustress* |
165 |
& (ustress(i,j,bi,bj)+ustress(i-1,j,bi,bj)) |
166 |
& *exf_half*maskW(i,j,ks,bi,bj) |
167 |
enddo |
168 |
enddo |
169 |
ENDIF |
170 |
|
171 |
#ifdef ALLOW_AUTODIFF_TAMC |
172 |
CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
173 |
#endif |
174 |
do j = jmin,jmax |
175 |
do i = imin,imax |
176 |
c Meridional wind stress. |
177 |
if (vstress(i,j,bi,bj).gt.windstressmax) then |
178 |
vstress(i,j,bi,bj)=windstressmax |
179 |
endif |
180 |
enddo |
181 |
enddo |
182 |
#ifdef ALLOW_AUTODIFF_TAMC |
183 |
CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
184 |
#endif |
185 |
do j = jmin,jmax |
186 |
do i = imin,imax |
187 |
if (vstress(i,j,bi,bj).lt.-windstressmax) then |
188 |
vstress(i,j,bi,bj)=-windstressmax |
189 |
endif |
190 |
enddo |
191 |
enddo |
192 |
IF ( stressIsOnCgrid ) THEN |
193 |
do j = jmin+1,jmax |
194 |
do i = imin,imax |
195 |
fv(i,j,bi,bj) = exf_outscal_vstress*vstress(i,j,bi,bj) |
196 |
enddo |
197 |
enddo |
198 |
ELSE |
199 |
do j = jmin+1,jmax |
200 |
do i = imin,imax |
201 |
c Shift wind stresses calculated at C-points to W/S points |
202 |
fv(i,j,bi,bj) = exf_outscal_vstress* |
203 |
& (vstress(i,j,bi,bj)+vstress(i,j-1,bi,bj)) |
204 |
& *exf_half*maskS(i,j,ks,bi,bj) |
205 |
enddo |
206 |
enddo |
207 |
ENDIF |
208 |
|
209 |
#ifdef SHORTWAVE_HEATING |
210 |
c Short wave radiative flux. |
211 |
do j = jmin,jmax |
212 |
do i = imin,imax |
213 |
qsw(i,j,bi,bj) = exf_outscal_swflux*swflux(i,j,bi,bj) |
214 |
enddo |
215 |
enddo |
216 |
#endif |
217 |
|
218 |
#ifdef ALLOW_CLIMSST_RELAXATION |
219 |
do j = jmin,jmax |
220 |
do i = imin,imax |
221 |
sst(i,j,bi,bj) = exf_outscal_sst*climsst(i,j,bi,bj) |
222 |
enddo |
223 |
enddo |
224 |
#endif |
225 |
|
226 |
#ifdef ALLOW_CLIMSSS_RELAXATION |
227 |
do j = jmin,jmax |
228 |
do i = imin,imax |
229 |
sss(i,j,bi,bj) = exf_outscal_sss*climsss(i,j,bi,bj) |
230 |
enddo |
231 |
enddo |
232 |
#endif |
233 |
|
234 |
#ifdef ATMOSPHERIC_LOADING |
235 |
do j = jmin,jmax |
236 |
do i = imin,imax |
237 |
pload(i,j,bi,bj)=exf_outscal_apressure*apressure(i,j,bi,bj) |
238 |
enddo |
239 |
enddo |
240 |
#endif |
241 |
|
242 |
#ifdef EXF_ALLOW_SEAICE_RELAX |
243 |
do j = jmin,jmax |
244 |
do i = imin,imax |
245 |
obsSIce(i,j,bi,bj) = |
246 |
& exf_outscal_areamask*areamask(i,j,bi,bj) |
247 |
obsSIce(I,J,bi,bj) = |
248 |
& MIN(MAX(obsSIce(I,J,bi,bj), 0.d0 ), 1.d0) |
249 |
enddo |
250 |
enddo |
251 |
#endif |
252 |
|
253 |
ENDDO |
254 |
ENDDO |
255 |
|
256 |
c Update the tile edges. |
257 |
|
258 |
_EXCH_XY_RS( qnet, mythid ) |
259 |
_EXCH_XY_RS( empmr, mythid ) |
260 |
CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid) |
261 |
#ifdef SHORTWAVE_HEATING |
262 |
_EXCH_XY_RS( qsw, mythid ) |
263 |
#endif |
264 |
#ifdef ALLOW_CLIMSST_RELAXATION |
265 |
_EXCH_XY_RS( sst, mythid ) |
266 |
#endif |
267 |
#ifdef ALLOW_CLIMSSS_RELAXATION |
268 |
_EXCH_XY_RS( sss, mythid ) |
269 |
#endif |
270 |
#ifdef ATMOSPHERIC_LOADING |
271 |
_EXCH_XY_RS( pload, mythid ) |
272 |
#endif |
273 |
#ifdef EXF_ALLOW_SEAICE_RELAX |
274 |
_EXCH_XY_RS( obsSIce, mythid) |
275 |
#endif |
276 |
|
277 |
RETURN |
278 |
END |