/[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.18 - (show annotations) (download)
Tue Apr 17 18:53:32 2007 UTC (17 years, 2 months ago) by jmc
Branch: MAIN
Changes since 1.17: +5 -3 lines
try to fix wind-stress location: defined on A-grid (if ALLOW_ATM_WIND &
 ALLOW_BULKFORMULAE or USE_EXF_INTERPOLATION), otherwise, defined on C-grid

1 C $Header: /u/gcmpack/MITgcm/pkg/exf/exf_mapfields.F,v 1.17 2007/04/16 23:27:21 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 #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 (ALLOW_ATM_WIND) ) \
159 || defined (USE_EXF_INTERPOLATION) )
160 c Shift wind stresses calculated at C-points 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))/2.*
163 & maskW(i,j,1,bi,bj)
164 #else
165 fu(i,j,bi,bj) = exf_outscal_ustress*ustress(i,j,bi,bj)
166 #endif
167 enddo
168 enddo
169
170 #ifdef ALLOW_AUTODIFF_TAMC
171 CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
172 #endif
173 do j = jmin,jmax
174 do i = imin,imax
175 c Meridional wind stress.
176 if (vstress(i,j,bi,bj).gt.windstressmax) then
177 vstress(i,j,bi,bj)=windstressmax
178 endif
179 enddo
180 enddo
181 #ifdef ALLOW_AUTODIFF_TAMC
182 CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
183 #endif
184 do j = jmin,jmax
185 do i = imin,imax
186 if (vstress(i,j,bi,bj).lt.-windstressmax) then
187 vstress(i,j,bi,bj)=-windstressmax
188 endif
189 enddo
190 enddo
191 do j = jmin+1,jmax
192 do i = imin,imax
193 #if ( ( defined (ALLOW_BULKFORMULAE) && defined (ALLOW_ATM_WIND) ) \
194 || defined (USE_EXF_INTERPOLATION) )
195 c Shift wind stresses calculated at C-points to W/S points
196 fv(i,j,bi,bj) = exf_outscal_vstress*
197 & (vstress(i,j,bi,bj)+vstress(i,j-1,bi,bj))/2.*
198 & maskS(i,j,1,bi,bj)
199 #else
200 fv(i,j,bi,bj) = exf_outscal_vstress*vstress(i,j,bi,bj)
201 #endif
202 enddo
203 enddo
204
205 #ifdef SHORTWAVE_HEATING
206 c Short wave radiative flux.
207 do j = jmin,jmax
208 do i = imin,imax
209 qsw(i,j,bi,bj) = exf_outscal_swflux*swflux(i,j,bi,bj)
210 enddo
211 enddo
212 #endif
213
214 #ifdef ALLOW_CLIMSST_RELAXATION
215 do j = jmin,jmax
216 do i = imin,imax
217 sst(i,j,bi,bj) = exf_outscal_sst*climsst(i,j,bi,bj)
218 enddo
219 enddo
220 #endif
221
222 #ifdef ALLOW_CLIMSSS_RELAXATION
223 do j = jmin,jmax
224 do i = imin,imax
225 sss(i,j,bi,bj) = exf_outscal_sss*climsss(i,j,bi,bj)
226 enddo
227 enddo
228 #endif
229
230 #ifdef ATMOSPHERIC_LOADING
231 do j = jmin,jmax
232 do i = imin,imax
233 pload(i,j,bi,bj)=exf_outscal_apressure*apressure(i,j,bi,bj)
234 enddo
235 enddo
236 #endif
237
238 enddo
239 enddo
240
241 c Update the tile edges.
242
243 _EXCH_XY_R4( qnet, mythid )
244 _EXCH_XY_R4( empmr, mythid )
245 c _EXCH_XY_R4( fu, mythid )
246 c _EXCH_XY_R4( fv, mythid )
247 CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid)
248 #ifdef SHORTWAVE_HEATING
249 _EXCH_XY_R4( qsw, mythid )
250 #endif
251 #ifdef ALLOW_CLIMSST_RELAXATION
252 _EXCH_XY_R4( sst, mythid )
253 #endif
254 #ifdef ALLOW_CLIMSSS_RELAXATION
255 _EXCH_XY_R4( sss, mythid )
256 #endif
257 #ifdef ATMOSPHERIC_LOADING
258 _EXCH_XY_R4( pload, mythid )
259 #endif
260
261 end

  ViewVC Help
Powered by ViewVC 1.1.22