/[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.23 - (show annotations) (download)
Tue Apr 28 18:15:33 2009 UTC (15 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62c, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62w, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint63g, checkpoint62, checkpoint63, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint62b, checkpoint61n, checkpoint61q, checkpoint61o, checkpoint61m, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.22: +7 -7 lines
change macros (EXCH & GLOBAL_SUM/MAX) sufix _R4/_R8 to _RS/_RL
 when applied to _RS/_RL variable

1 C $Header: /u/gcmpack/MITgcm/pkg/exf/exf_mapfields.F,v 1.22 2007/10/01 13:36:05 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,k
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 DO bj = myByLo(myThid),myByHi(myThid)
82 DO bi = myBxLo(myThid),myBxHi(myThid)
83
84 #ifdef ALLOW_AUTODIFF_TAMC
85 act1 = bi - myBxLo(myThid)
86 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
87 act2 = bj - myByLo(myThid)
88 max2 = myByHi(myThid) - myByLo(myThid) + 1
89 act3 = myThid - 1
90 max3 = nTx*nTy
91 act4 = ikey_dynamics - 1
92 ikey = (act1 + 1) + act2*max1
93 & + act3*max1*max2
94 & + act4*max1*max2*max3
95 #endif /* ALLOW_AUTODIFF_TAMC */
96
97 c Heat flux.
98 do j = jmin,jmax
99 do i = imin,imax
100 qnet(i,j,bi,bj) = exf_outscal_hflux*hflux(i,j,bi,bj)
101 enddo
102 enddo
103 if ( hfluxfile .EQ. ' ' ) then
104 do j = jmin,jmax
105 do i = imin,imax
106 qnet(i,j,bi,bj) = qnet(i,j,bi,bj) -
107 & exf_outscal_hflux * ( hflux_exfremo_intercept +
108 & hflux_exfremo_slope*(mytime-starttime) )
109 enddo
110 enddo
111 endif
112
113 c Salt flux.
114 do j = jmin,jmax
115 do i = imin,imax
116 EmPmR(i,j,bi,bj)= exf_outscal_sflux*sflux(i,j,bi,bj)
117 & *rhoConstFresh
118 enddo
119 enddo
120 if ( sfluxfile .EQ. ' ' ) then
121 do j = jmin,jmax
122 do i = imin,imax
123 EmPmR(i,j,bi,bj) = EmPmR(i,j,bi,bj) - rhoConstFresh*
124 & exf_outscal_sflux * ( sflux_exfremo_intercept +
125 & sflux_exfremo_slope*(mytime-starttime) )
126 enddo
127 enddo
128 endif
129
130 #ifdef ALLOW_AUTODIFF_TAMC
131 CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
132 #endif
133 do j = jmin,jmax
134 do i = imin,imax
135 c Zonal wind stress.
136 if (ustress(i,j,bi,bj).gt.windstressmax) then
137 ustress(i,j,bi,bj)=windstressmax
138 endif
139 enddo
140 enddo
141 #ifdef ALLOW_AUTODIFF_TAMC
142 CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
143 #endif
144 do j = jmin,jmax
145 do i = imin,imax
146 if (ustress(i,j,bi,bj).lt.-windstressmax) then
147 ustress(i,j,bi,bj)=-windstressmax
148 endif
149 enddo
150 enddo
151 IF ( stressIsOnCgrid ) THEN
152 do j = jmin,jmax
153 do i = imin+1,imax
154 fu(i,j,bi,bj) = exf_outscal_ustress*ustress(i,j,bi,bj)
155 enddo
156 enddo
157 ELSE
158 do j = jmin,jmax
159 do i = imin+1,imax
160 c Shift wind stresses calculated at Grid-center to W/S points
161 fu(i,j,bi,bj) = exf_outscal_ustress*
162 & (ustress(i,j,bi,bj)+ustress(i-1,j,bi,bj))
163 & *exf_half*maskW(i,j,1,bi,bj)
164 enddo
165 enddo
166 ENDIF
167
168 #ifdef ALLOW_AUTODIFF_TAMC
169 CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
170 #endif
171 do j = jmin,jmax
172 do i = imin,imax
173 c Meridional wind stress.
174 if (vstress(i,j,bi,bj).gt.windstressmax) then
175 vstress(i,j,bi,bj)=windstressmax
176 endif
177 enddo
178 enddo
179 #ifdef ALLOW_AUTODIFF_TAMC
180 CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
181 #endif
182 do j = jmin,jmax
183 do i = imin,imax
184 if (vstress(i,j,bi,bj).lt.-windstressmax) then
185 vstress(i,j,bi,bj)=-windstressmax
186 endif
187 enddo
188 enddo
189 IF ( stressIsOnCgrid ) THEN
190 do j = jmin+1,jmax
191 do i = imin,imax
192 fv(i,j,bi,bj) = exf_outscal_vstress*vstress(i,j,bi,bj)
193 enddo
194 enddo
195 ELSE
196 do j = jmin+1,jmax
197 do i = imin,imax
198 c Shift wind stresses calculated at C-points to W/S points
199 fv(i,j,bi,bj) = exf_outscal_vstress*
200 & (vstress(i,j,bi,bj)+vstress(i,j-1,bi,bj))
201 & *exf_half*maskS(i,j,1,bi,bj)
202 enddo
203 enddo
204 ENDIF
205
206 #ifdef SHORTWAVE_HEATING
207 c Short wave radiative flux.
208 do j = jmin,jmax
209 do i = imin,imax
210 qsw(i,j,bi,bj) = exf_outscal_swflux*swflux(i,j,bi,bj)
211 enddo
212 enddo
213 #endif
214
215 #ifdef ALLOW_CLIMSST_RELAXATION
216 do j = jmin,jmax
217 do i = imin,imax
218 sst(i,j,bi,bj) = exf_outscal_sst*climsst(i,j,bi,bj)
219 enddo
220 enddo
221 #endif
222
223 #ifdef ALLOW_CLIMSSS_RELAXATION
224 do j = jmin,jmax
225 do i = imin,imax
226 sss(i,j,bi,bj) = exf_outscal_sss*climsss(i,j,bi,bj)
227 enddo
228 enddo
229 #endif
230
231 #ifdef ATMOSPHERIC_LOADING
232 do j = jmin,jmax
233 do i = imin,imax
234 pload(i,j,bi,bj)=exf_outscal_apressure*apressure(i,j,bi,bj)
235 enddo
236 enddo
237 #endif
238
239 ENDDO
240 ENDDO
241
242 c Update the tile edges.
243
244 _EXCH_XY_RS( qnet, mythid )
245 _EXCH_XY_RS( empmr, mythid )
246 CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid)
247 #ifdef SHORTWAVE_HEATING
248 _EXCH_XY_RS( qsw, mythid )
249 #endif
250 #ifdef ALLOW_CLIMSST_RELAXATION
251 _EXCH_XY_RS( sst, mythid )
252 #endif
253 #ifdef ALLOW_CLIMSSS_RELAXATION
254 _EXCH_XY_RS( sss, mythid )
255 #endif
256 #ifdef ATMOSPHERIC_LOADING
257 _EXCH_XY_RS( pload, mythid )
258 #endif
259
260 RETURN
261 END

  ViewVC Help
Powered by ViewVC 1.1.22