1 |
C $Header: $ |
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 |
#include "EXF_CLIM_PARAM.h" |
59 |
#include "EXF_CLIM_FIELDS.h" |
60 |
#ifdef ALLOW_AUTODIFF_TAMC |
61 |
# include "tamc.h" |
62 |
# include "tamc_keys.h" |
63 |
#endif |
64 |
c == routine arguments == |
65 |
|
66 |
c mythid - thread number for this instance of the routine. |
67 |
|
68 |
integer mythid |
69 |
integer myiter |
70 |
_RL mytime |
71 |
|
72 |
c == local variables == |
73 |
|
74 |
integer bi,bj |
75 |
integer i,j,k |
76 |
integer jtlo |
77 |
integer jthi |
78 |
integer itlo |
79 |
integer ithi |
80 |
integer jmin |
81 |
integer jmax |
82 |
integer imin |
83 |
integer imax |
84 |
|
85 |
c == end of interface == |
86 |
|
87 |
jtlo = mybylo(mythid) |
88 |
jthi = mybyhi(mythid) |
89 |
itlo = mybxlo(mythid) |
90 |
ithi = mybxhi(mythid) |
91 |
jmin = 1-oly |
92 |
jmax = sny+oly |
93 |
imin = 1-olx |
94 |
imax = snx+olx |
95 |
|
96 |
do bj = jtlo,jthi |
97 |
do bi = itlo,ithi |
98 |
|
99 |
#ifdef ALLOW_AUTODIFF_TAMC |
100 |
act1 = bi - myBxLo(myThid) |
101 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
102 |
act2 = bj - myByLo(myThid) |
103 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
104 |
act3 = myThid - 1 |
105 |
max3 = nTx*nTy |
106 |
act4 = ikey_dynamics - 1 |
107 |
ikey = (act1 + 1) + act2*max1 |
108 |
& + act3*max1*max2 |
109 |
& + act4*max1*max2*max3 |
110 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
111 |
|
112 |
do j = jmin,jmax |
113 |
do i = imin,imax |
114 |
c Heat flux. |
115 |
qnet(i,j,bi,bj) = exf_outscal_hflux*hflux(i,j,bi,bj) |
116 |
if ( hfluxfile .EQ. ' ' ) |
117 |
& qnet(i,j,bi,bj) = qnet(i,j,bi,bj) - |
118 |
& exf_outscal_hflux * ( hflux_exfremo_intercept + |
119 |
& hflux_exfremo_slope*(mytime-starttime) ) |
120 |
enddo |
121 |
enddo |
122 |
|
123 |
|
124 |
do j = jmin,jmax |
125 |
do i = imin,imax |
126 |
c Salt flux. |
127 |
empmr(i,j,bi,bj)= exf_outscal_sflux*sflux(i,j,bi,bj) |
128 |
if ( sfluxfile .EQ. ' ' ) |
129 |
& empmr(i,j,bi,bj) = empmr(i,j,bi,bj) - |
130 |
& exf_outscal_sflux * ( sflux_exfremo_intercept + |
131 |
& sflux_exfremo_slope*(mytime-starttime) ) |
132 |
enddo |
133 |
enddo |
134 |
|
135 |
#ifdef ALLOW_AUTODIFF_TAMC |
136 |
CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
137 |
#endif |
138 |
do j = jmin,jmax |
139 |
do i = imin,imax |
140 |
c Zonal wind stress. |
141 |
if (ustress(i,j,bi,bj).gt.windstressmax) then |
142 |
ustress(i,j,bi,bj)=windstressmax |
143 |
endif |
144 |
enddo |
145 |
enddo |
146 |
#ifdef ALLOW_AUTODIFF_TAMC |
147 |
CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
148 |
#endif |
149 |
do j = jmin,jmax |
150 |
do i = imin,imax |
151 |
if (ustress(i,j,bi,bj).lt.-windstressmax) then |
152 |
ustress(i,j,bi,bj)=-windstressmax |
153 |
endif |
154 |
enddo |
155 |
enddo |
156 |
do j = jmin,jmax |
157 |
do i = imin+1,imax |
158 |
#if (defined (ALLOW_BULKFORMULAE) || defined (USE_EXF_INTERPOLATION)) |
159 |
c Shift wind stresses calculated at C-points to W/S points |
160 |
fu(i,j,bi,bj) = exf_outscal_ustress* |
161 |
& (ustress(i,j,bi,bj)+ustress(i-1,j,bi,bj))/2.* |
162 |
& maskW(i,j,1,bi,bj) |
163 |
#else |
164 |
fu(i,j,bi,bj) = exf_outscal_ustress*ustress(i,j,bi,bj) |
165 |
#endif |
166 |
enddo |
167 |
enddo |
168 |
|
169 |
#ifdef ALLOW_AUTODIFF_TAMC |
170 |
CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
171 |
#endif |
172 |
do j = jmin,jmax |
173 |
do i = imin,imax |
174 |
c Meridional wind stress. |
175 |
if (vstress(i,j,bi,bj).gt.windstressmax) then |
176 |
vstress(i,j,bi,bj)=windstressmax |
177 |
endif |
178 |
enddo |
179 |
enddo |
180 |
#ifdef ALLOW_AUTODIFF_TAMC |
181 |
CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
182 |
#endif |
183 |
do j = jmin,jmax |
184 |
do i = imin,imax |
185 |
if (vstress(i,j,bi,bj).lt.-windstressmax) then |
186 |
vstress(i,j,bi,bj)=-windstressmax |
187 |
endif |
188 |
enddo |
189 |
enddo |
190 |
do j = jmin+1,jmax |
191 |
do i = imin,imax |
192 |
#if (defined (ALLOW_BULKFORMULAE) || defined (USE_EXF_INTERPOLATION)) |
193 |
c Shift wind stresses calculated at C-points to W/S points |
194 |
fv(i,j,bi,bj) = exf_outscal_vstress* |
195 |
& (vstress(i,j,bi,bj)+vstress(i,j-1,bi,bj))/2.* |
196 |
& maskS(i,j,1,bi,bj) |
197 |
#else |
198 |
fv(i,j,bi,bj) = exf_outscal_vstress*vstress(i,j,bi,bj) |
199 |
#endif |
200 |
enddo |
201 |
enddo |
202 |
|
203 |
#ifdef SHORTWAVE_HEATING |
204 |
c Short wave radiative flux. |
205 |
do j = jmin,jmax |
206 |
do i = imin,imax |
207 |
qsw(i,j,bi,bj) = exf_outscal_swflux*swflux(i,j,bi,bj) |
208 |
enddo |
209 |
enddo |
210 |
#endif |
211 |
|
212 |
#ifdef ALLOW_CLIMSST_RELAXATION |
213 |
do j = jmin,jmax |
214 |
do i = imin,imax |
215 |
sst(i,j,bi,bj) = exf_outscal_sst*climsst(i,j,bi,bj) |
216 |
enddo |
217 |
enddo |
218 |
#endif |
219 |
|
220 |
#ifdef ALLOW_CLIMSSS_RELAXATION |
221 |
do j = jmin,jmax |
222 |
do i = imin,imax |
223 |
sss(i,j,bi,bj) = exf_outscal_sss*climsss(i,j,bi,bj) |
224 |
enddo |
225 |
enddo |
226 |
#endif |
227 |
|
228 |
#ifdef ATMOSPHERIC_LOADING |
229 |
do j = jmin,jmax |
230 |
do i = imin,imax |
231 |
pload(i,j,bi,bj)=exf_outscal_apressure*apressure(i,j,bi,bj) |
232 |
enddo |
233 |
enddo |
234 |
#endif |
235 |
|
236 |
enddo |
237 |
enddo |
238 |
|
239 |
c Update the tile edges. |
240 |
|
241 |
_EXCH_XY_R4( qnet, mythid ) |
242 |
_EXCH_XY_R4( empmr, mythid ) |
243 |
c _EXCH_XY_R4( fu, mythid ) |
244 |
c _EXCH_XY_R4( fv, mythid ) |
245 |
CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid) |
246 |
#ifdef SHORTWAVE_HEATING |
247 |
_EXCH_XY_R4( qsw, mythid ) |
248 |
#endif |
249 |
#ifdef ALLOW_CLIMSST_RELAXATION |
250 |
_EXCH_XY_R4( sst, mythid ) |
251 |
#endif |
252 |
#ifdef ALLOW_CLIMSSS_RELAXATION |
253 |
_EXCH_XY_R4( sss, mythid ) |
254 |
#endif |
255 |
#ifdef ATMOSPHERIC_LOADING |
256 |
_EXCH_XY_R4( pload, mythid ) |
257 |
#endif |
258 |
|
259 |
end |