/[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.14 - (show annotations) (download)
Mon Jan 2 21:17:02 2006 UTC (18 years, 5 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58
Changes since 1.13: +1 -25 lines
o Fix I/O inconsistency in pkg/rbcs: replace precFloat32 by readBinaryPrec
o Remove 3-dim. relaxation code from pkg/exf (now use only pkg/rbcs)
o Thanks to Tom Haine for testing!

1 c $Header: /u/gcmpack/MITgcm/pkg/exf/exf_mapfields.F,v 1.13 2005/12/13 19:46:46 heimbach Exp $
2
3 #include "EXF_OPTIONS.h"
4
5 subroutine exf_mapfields( mythid )
6
7 c ==================================================================
8 c SUBROUTINE exf_mapfields
9 c ==================================================================
10 c
11 c o Map external forcing fields (ustress, vstress, hflux, sflux,
12 c swflux, apressure, climsss, climsst, etc.) onto ocean model
13 c arrays (fu, fv, Qnet, EmPmR, Qsw, pload, sss, sst, etc.).
14 c This routine is included to separate the ocean state estimation
15 c tool as much as possible from the ocean model. Unit and sign
16 c conventions can be customized using variables exf_outscal_*,
17 c which are set in exf_readparms.F. See the header files
18 c exf_fields.h and FFIELDS.h for definitions of the various input
19 c and output fields and for default unit and sign convetions.
20 c
21 c started: Christian Eckert eckert@mit.edu 09-Aug-1999
22 c
23 c changed: Christian Eckert eckert@mit.edu 11-Jan-2000
24 c - Restructured the code in order to create a package
25 c for the MITgcmUV.
26 c
27 c Christian Eckert eckert@mit.edu 12-Feb-2000
28 c - Changed Routine names (package prefix: exf_)
29 c
30 c Patrick Heimbach, heimbach@mit.edu 06-May-2000
31 c - added and changed CPP flag structure for
32 c ALLOW_BULKFORMULAE, ALLOW_ATM_TEMP
33 c
34 c Patrick Heimbach, heimbach@mit.edu 23-May-2000
35 c - sign change of ustress/vstress incorporated into
36 c scaling factors exf_outscal_ust, exf_outscal_vst
37 c
38 c mods for pkg/seaice: menemenlis@jpl.nasa.gov 20-Dec-2002
39 c
40 c ==================================================================
41 c SUBROUTINE exf_mapfields
42 c ==================================================================
43
44 implicit none
45
46 c == global variables ==
47
48 #include "EEPARAMS.h"
49 #include "SIZE.h"
50 #include "FFIELDS.h"
51 #include "GRID.h"
52
53 #include "exf_param.h"
54 #include "exf_constants.h"
55 #include "exf_fields.h"
56 #include "exf_clim_param.h"
57 #include "exf_clim_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 integer mythid
67
68 c == local variables ==
69
70 integer bi,bj
71 integer i,j,k
72 integer jtlo
73 integer jthi
74 integer itlo
75 integer ithi
76 integer jmin
77 integer jmax
78 integer imin
79 integer imax
80
81 c == end of interface ==
82
83 jtlo = mybylo(mythid)
84 jthi = mybyhi(mythid)
85 itlo = mybxlo(mythid)
86 ithi = mybxhi(mythid)
87 jmin = 1-oly
88 jmax = sny+oly
89 imin = 1-olx
90 imax = snx+olx
91
92 do bj = jtlo,jthi
93 do bi = itlo,ithi
94
95 #ifdef ALLOW_AUTODIFF_TAMC
96 act1 = bi - myBxLo(myThid)
97 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
98 act2 = bj - myByLo(myThid)
99 max2 = myByHi(myThid) - myByLo(myThid) + 1
100 act3 = myThid - 1
101 max3 = nTx*nTy
102 act4 = ikey_dynamics - 1
103 ikey = (act1 + 1) + act2*max1
104 & + act3*max1*max2
105 & + act4*max1*max2*max3
106 #endif /* ALLOW_AUTODIFF_TAMC */
107
108 do j = jmin,jmax
109 do i = imin,imax
110 c Heat flux.
111 qnet(i,j,bi,bj) = exf_outscal_hflux*hflux(i,j,bi,bj)
112 enddo
113 enddo
114
115
116 do j = jmin,jmax
117 do i = imin,imax
118 c Salt flux.
119 empmr(i,j,bi,bj)= exf_outscal_sflux*sflux(i,j,bi,bj)
120 enddo
121 enddo
122
123 #ifdef ALLOW_AUTODIFF_TAMC
124 CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
125 #endif
126 do j = jmin,jmax
127 do i = imin,imax
128 c Zonal wind stress.
129 if (ustress(i,j,bi,bj).gt.windstressmax) then
130 ustress(i,j,bi,bj)=windstressmax
131 endif
132 enddo
133 enddo
134 #ifdef ALLOW_AUTODIFF_TAMC
135 CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
136 #endif
137 do j = jmin,jmax
138 do i = imin,imax
139 if (ustress(i,j,bi,bj).lt.-windstressmax) then
140 ustress(i,j,bi,bj)=-windstressmax
141 endif
142 enddo
143 enddo
144 do j = jmin,jmax
145 do i = imin+1,imax
146 #if (defined (ALLOW_BULKFORMULAE) || defined (USE_EXF_INTERPOLATION))
147 c Shift wind stresses calculated at C-points to W/S points
148 fu(i,j,bi,bj) = exf_outscal_ustress*
149 & (ustress(i,j,bi,bj)+ustress(i-1,j,bi,bj))/2.*
150 & maskW(i,j,1,bi,bj)
151 #else
152 fu(i,j,bi,bj) = exf_outscal_ustress*ustress(i,j,bi,bj)
153 #endif
154 enddo
155 enddo
156
157 #ifdef ALLOW_AUTODIFF_TAMC
158 CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
159 #endif
160 do j = jmin,jmax
161 do i = imin,imax
162 c Meridional wind stress.
163 if (vstress(i,j,bi,bj).gt.windstressmax) then
164 vstress(i,j,bi,bj)=windstressmax
165 endif
166 enddo
167 enddo
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 if (vstress(i,j,bi,bj).lt.-windstressmax) then
174 vstress(i,j,bi,bj)=-windstressmax
175 endif
176 enddo
177 enddo
178 do j = jmin+1,jmax
179 do i = imin,imax
180 #if (defined (ALLOW_BULKFORMULAE) || defined (USE_EXF_INTERPOLATION))
181 c Shift wind stresses calculated at C-points to W/S points
182 fv(i,j,bi,bj) = exf_outscal_vstress*
183 & (vstress(i,j,bi,bj)+vstress(i,j-1,bi,bj))/2.*
184 & maskS(i,j,1,bi,bj)
185 #else
186 fv(i,j,bi,bj) = exf_outscal_vstress*vstress(i,j,bi,bj)
187 #endif
188 enddo
189 enddo
190
191 #ifdef SHORTWAVE_HEATING
192 c Short wave radiative flux.
193 do j = jmin,jmax
194 do i = imin,imax
195 qsw(i,j,bi,bj) = exf_outscal_swflux*swflux(i,j,bi,bj)
196 enddo
197 enddo
198 #endif
199
200 #ifdef ALLOW_CLIMSST_RELAXATION
201 do j = jmin,jmax
202 do i = imin,imax
203 sst(i,j,bi,bj) = exf_outscal_sst*climsst(i,j,bi,bj)
204 enddo
205 enddo
206 #endif
207
208 #ifdef ALLOW_CLIMSSS_RELAXATION
209 do j = jmin,jmax
210 do i = imin,imax
211 sss(i,j,bi,bj) = exf_outscal_sss*climsss(i,j,bi,bj)
212 enddo
213 enddo
214 #endif
215
216 #ifdef ATMOSPHERIC_LOADING
217 do j = jmin,jmax
218 do i = imin,imax
219 pload(i,j,bi,bj)=exf_outscal_apressure*apressure(i,j,bi,bj)
220 enddo
221 enddo
222 #endif
223
224 enddo
225 enddo
226
227 c Update the tile edges.
228
229 _EXCH_XY_R4( qnet, mythid )
230 _EXCH_XY_R4( empmr, mythid )
231 c _EXCH_XY_R4( fu, mythid )
232 c _EXCH_XY_R4( fv, mythid )
233 CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid)
234 #ifdef SHORTWAVE_HEATING
235 _EXCH_XY_R4( qsw, mythid )
236 #endif
237 #ifdef ALLOW_CLIMSST_RELAXATION
238 _EXCH_XY_R4( sst, mythid )
239 #endif
240 #ifdef ALLOW_CLIMSSS_RELAXATION
241 _EXCH_XY_R4( sss, mythid )
242 #endif
243 #ifdef ALLOW_CLIMTEMP_RELAXATION
244 _EXCH_XYZ_R4( thetaStar, mythid )
245 #endif
246 #ifdef ALLOW_CLIMSALT_RELAXATION
247 _EXCH_XYZ_R4( saltStar, mythid )
248 #endif
249 #ifdef ATMOSPHERIC_LOADING
250 _EXCH_XY_R4( pload, mythid )
251 #endif
252
253 end

  ViewVC Help
Powered by ViewVC 1.1.22