1 |
c $Header: /u/u3/gcmpack/MITgcm/pkg/exf/exf_mapfields.F,v 1.6.2.1 2003/10/02 18:30:07 adcroft Exp $ |
2 |
|
3 |
#include "EXF_OPTIONS.h" |
4 |
|
5 |
|
6 |
subroutine exf_mapfields( 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 "FFIELDS.h" |
52 |
#include "exf_param.h" |
53 |
#include "exf_constants.h" |
54 |
#include "exf_fields.h" |
55 |
#include "exf_clim_fields.h" |
56 |
#ifdef ALLOW_AUTODIFF_TAMC |
57 |
# include "tamc.h" |
58 |
# include "tamc_keys.h" |
59 |
#endif |
60 |
c == routine arguments == |
61 |
|
62 |
c mythid - thread number for this instance of the routine. |
63 |
|
64 |
integer mythid |
65 |
|
66 |
c == local variables == |
67 |
|
68 |
integer bi,bj |
69 |
integer i,j |
70 |
integer jtlo |
71 |
integer jthi |
72 |
integer itlo |
73 |
integer ithi |
74 |
integer jmin |
75 |
integer jmax |
76 |
integer imin |
77 |
integer imax |
78 |
|
79 |
c == end of interface == |
80 |
|
81 |
jtlo = mybylo(mythid) |
82 |
jthi = mybyhi(mythid) |
83 |
itlo = mybxlo(mythid) |
84 |
ithi = mybxhi(mythid) |
85 |
jmin = 1-oly |
86 |
jmax = sny+oly |
87 |
imin = 1-olx |
88 |
imax = snx+olx |
89 |
|
90 |
do bj = jtlo,jthi |
91 |
do bi = itlo,ithi |
92 |
|
93 |
#ifdef ALLOW_AUTODIFF_TAMC |
94 |
act1 = bi - myBxLo(myThid) |
95 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
96 |
act2 = bj - myByLo(myThid) |
97 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
98 |
act3 = myThid - 1 |
99 |
max3 = nTx*nTy |
100 |
act4 = ikey_dynamics - 1 |
101 |
ikey = (act1 + 1) + act2*max1 |
102 |
& + act3*max1*max2 |
103 |
& + act4*max1*max2*max3 |
104 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
105 |
|
106 |
do j = jmin,jmax |
107 |
do i = imin,imax |
108 |
c Heat flux. |
109 |
qnet(i,j,bi,bj) = exf_outscal_hflux*hflux(i,j,bi,bj) |
110 |
enddo |
111 |
enddo |
112 |
|
113 |
|
114 |
do j = jmin,jmax |
115 |
do i = imin,imax |
116 |
c Salt flux. |
117 |
empmr(i,j,bi,bj)= exf_outscal_sflux*sflux(i,j,bi,bj) |
118 |
enddo |
119 |
enddo |
120 |
|
121 |
#ifdef ALLOW_AUTODIFF_TAMC |
122 |
CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
123 |
#endif |
124 |
do j = jmin,jmax |
125 |
do i = imin,imax |
126 |
c Zonal wind stress. |
127 |
if (ustress(i,j,bi,bj).gt.2.0D0) then |
128 |
ustress(i,j,bi,bj)=2.0D0 |
129 |
endif |
130 |
enddo |
131 |
enddo |
132 |
#ifdef ALLOW_AUTODIFF_TAMC |
133 |
CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
134 |
#endif |
135 |
do j = jmin,jmax |
136 |
do i = imin,imax |
137 |
if (ustress(i,j,bi,bj).lt.-2.0D0) then |
138 |
ustress(i,j,bi,bj)=-2.0D0 |
139 |
endif |
140 |
enddo |
141 |
enddo |
142 |
do j = jmin,jmax |
143 |
do i = imin,imax |
144 |
fu(i,j,bi,bj) = exf_outscal_ustress*ustress(i,j,bi,bj) |
145 |
enddo |
146 |
enddo |
147 |
|
148 |
#ifdef ALLOW_AUTODIFF_TAMC |
149 |
CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
150 |
#endif |
151 |
do j = jmin,jmax |
152 |
do i = imin,imax |
153 |
c Meridional wind stress. |
154 |
if (vstress(i,j,bi,bj).gt.2.0D0) then |
155 |
vstress(i,j,bi,bj)=2.0D0 |
156 |
endif |
157 |
enddo |
158 |
enddo |
159 |
#ifdef ALLOW_AUTODIFF_TAMC |
160 |
CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
161 |
#endif |
162 |
do j = jmin,jmax |
163 |
do i = imin,imax |
164 |
if (vstress(i,j,bi,bj).lt.-2.0D0) then |
165 |
vstress(i,j,bi,bj)=-2.0D0 |
166 |
endif |
167 |
enddo |
168 |
enddo |
169 |
do j = jmin,jmax |
170 |
do i = imin,imax |
171 |
fv(i,j,bi,bj) = exf_outscal_vstress*vstress(i,j,bi,bj) |
172 |
enddo |
173 |
enddo |
174 |
|
175 |
#ifdef SHORTWAVE_HEATING |
176 |
c Short wave radiative flux. |
177 |
do j = jmin,jmax |
178 |
do i = imin,imax |
179 |
qsw(i,j,bi,bj) = exf_outscal_swflux*swflux(i,j,bi,bj) |
180 |
enddo |
181 |
enddo |
182 |
#endif |
183 |
|
184 |
#ifdef ALLOW_CLIMSST_RELAXATION |
185 |
do j = jmin,jmax |
186 |
do i = imin,imax |
187 |
sst(i,j,bi,bj) = exf_outscal_sst*climsst(i,j,bi,bj) |
188 |
enddo |
189 |
enddo |
190 |
#endif |
191 |
|
192 |
#ifdef ALLOW_CLIMSSS_RELAXATION |
193 |
do j = jmin,jmax |
194 |
do i = imin,imax |
195 |
sss(i,j,bi,bj) = exf_outscal_sss*climsss(i,j,bi,bj) |
196 |
enddo |
197 |
enddo |
198 |
#endif |
199 |
|
200 |
#ifdef ATMOSPHERIC_LOADING |
201 |
do j = jmin,jmax |
202 |
do i = imin,imax |
203 |
pload(i,j,bi,bj)=exf_outscal_apressure*apressure(i,j,bi,bj) |
204 |
enddo |
205 |
enddo |
206 |
#endif |
207 |
|
208 |
enddo |
209 |
enddo |
210 |
|
211 |
c Update the tile edges. |
212 |
|
213 |
_EXCH_XY_R4( qnet, mythid ) |
214 |
_EXCH_XY_R4( empmr, mythid ) |
215 |
c _EXCH_XY_R4( fu, mythid ) |
216 |
c _EXCH_XY_R4( fv, mythid ) |
217 |
CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid) |
218 |
#ifdef SHORTWAVE_HEATING |
219 |
_EXCH_XY_R4( qsw, mythid ) |
220 |
#endif |
221 |
#ifdef ALLOW_CLIMSST_RELAXATION |
222 |
_EXCH_XY_R4( sst, mythid ) |
223 |
#endif |
224 |
#ifdef ALLOW_CLIMSSS_RELAXATION |
225 |
_EXCH_XY_R4( sss, mythid ) |
226 |
#endif |
227 |
#ifdef ATMOSPHERIC_LOADING |
228 |
_EXCH_XY_R4( pload, mythid ) |
229 |
#endif |
230 |
|
231 |
end |