/[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.20 - (show annotations) (download)
Mon May 14 19:34:57 2007 UTC (17 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59c, checkpoint59b
Changes since 1.19: +40 -48 lines
- implement A-grid / C-grid selection for wind-stress input files ;
- call the appropriate EXCH ; add some consistency check .

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

  ViewVC Help
Powered by ViewVC 1.1.22