/[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.26 - (show annotations) (download)
Thu Jan 10 22:59:24 2013 UTC (11 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64c, checkpoint64e, checkpoint64d, checkpoint64f
Changes since 1.25: +165 -160 lines
for diagnostics purpose, fill up Qsw array even when SHORTWAVE_HEATING is #undef

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

  ViewVC Help
Powered by ViewVC 1.1.22