/[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.11 - (show annotations) (download)
Fri Jul 16 01:20:57 2004 UTC (19 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57g_pre, checkpoint57b_post, checkpoint57g_post, checkpoint56b_post, checkpoint54d_post, checkpoint54e_post, checkpoint57d_post, checkpoint57i_post, checkpoint55, checkpoint57, checkpoint56, checkpoint54f_post, checkpoint55i_post, checkpoint55c_post, checkpoint57f_post, checkpoint57a_post, checkpoint57h_pre, checkpoint57h_post, checkpoint55g_post, checkpoint57c_post, checkpoint55d_post, checkpoint55d_pre, checkpoint57c_pre, checkpoint55j_post, checkpoint55h_post, checkpoint57e_post, checkpoint55b_post, checkpoint55f_post, eckpoint57e_pre, checkpoint56a_post, checkpoint57h_done, checkpoint57f_pre, checkpoint56c_post, checkpoint57a_pre, checkpoint55a_post, checkpoint55e_post, checkpoint54c_post
Changes since 1.10: +7 -1 lines
change Qnet to always be the net heat flux, (+upward).

1 c $Header: /u/gcmpack/MITgcm/pkg/exf/exf_mapfields.F,v 1.10 2004/04/19 22:30:46 mlosch 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_fields.h"
57 #ifdef ALLOW_AUTODIFF_TAMC
58 # include "tamc.h"
59 # include "tamc_keys.h"
60 #endif
61 c == routine arguments ==
62
63 c mythid - thread number for this instance of the routine.
64
65 integer mythid
66
67 c == local variables ==
68
69 integer bi,bj
70 integer i,j
71 integer jtlo
72 integer jthi
73 integer itlo
74 integer ithi
75 integer jmin
76 integer jmax
77 integer imin
78 integer imax
79
80 c == end of interface ==
81
82 jtlo = mybylo(mythid)
83 jthi = mybyhi(mythid)
84 itlo = mybxlo(mythid)
85 ithi = mybxhi(mythid)
86 jmin = 1-oly
87 jmax = sny+oly
88 imin = 1-olx
89 imax = snx+olx
90
91 do bj = jtlo,jthi
92 do bi = itlo,ithi
93
94 #ifdef ALLOW_AUTODIFF_TAMC
95 act1 = bi - myBxLo(myThid)
96 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
97 act2 = bj - myByLo(myThid)
98 max2 = myByHi(myThid) - myByLo(myThid) + 1
99 act3 = myThid - 1
100 max3 = nTx*nTy
101 act4 = ikey_dynamics - 1
102 ikey = (act1 + 1) + act2*max1
103 & + act3*max1*max2
104 & + act4*max1*max2*max3
105 #endif /* ALLOW_AUTODIFF_TAMC */
106
107 do j = jmin,jmax
108 do i = imin,imax
109 c Heat flux.
110 qnet(i,j,bi,bj) = exf_outscal_hflux*hflux(i,j,bi,bj)
111 enddo
112 enddo
113
114
115 do j = jmin,jmax
116 do i = imin,imax
117 c Salt flux.
118 empmr(i,j,bi,bj)= exf_outscal_sflux*sflux(i,j,bi,bj)
119 enddo
120 enddo
121
122 #ifdef ALLOW_AUTODIFF_TAMC
123 CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
124 #endif
125 do j = jmin,jmax
126 do i = imin,imax
127 c Zonal wind stress.
128 if (ustress(i,j,bi,bj).gt.windstressmax) then
129 ustress(i,j,bi,bj)=windstressmax
130 endif
131 enddo
132 enddo
133 #ifdef ALLOW_AUTODIFF_TAMC
134 CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
135 #endif
136 do j = jmin,jmax
137 do i = imin,imax
138 if (ustress(i,j,bi,bj).lt.-windstressmax) then
139 ustress(i,j,bi,bj)=-windstressmax
140 endif
141 enddo
142 enddo
143 do j = jmin,jmax
144 do i = imin+1,imax
145 #if (defined (ALLOW_BULKFORMULAE) || defined (USE_EXF_INTERPOLATION))
146 c Shift wind stresses calculated at C-points to W/S points
147 fu(i,j,bi,bj) = exf_outscal_ustress*
148 & (ustress(i,j,bi,bj)+ustress(i-1,j,bi,bj))/2.*
149 & maskW(i,j,1,bi,bj)
150 #else
151 fu(i,j,bi,bj) = exf_outscal_ustress*ustress(i,j,bi,bj)
152 #endif
153 enddo
154 enddo
155
156 #ifdef ALLOW_AUTODIFF_TAMC
157 CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
158 #endif
159 do j = jmin,jmax
160 do i = imin,imax
161 c Meridional wind stress.
162 if (vstress(i,j,bi,bj).gt.windstressmax) then
163 vstress(i,j,bi,bj)=windstressmax
164 endif
165 enddo
166 enddo
167 #ifdef ALLOW_AUTODIFF_TAMC
168 CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
169 #endif
170 do j = jmin,jmax
171 do i = imin,imax
172 if (vstress(i,j,bi,bj).lt.-windstressmax) then
173 vstress(i,j,bi,bj)=-windstressmax
174 endif
175 enddo
176 enddo
177 do j = jmin+1,jmax
178 do i = imin,imax
179 #if (defined (ALLOW_BULKFORMULAE) || defined (USE_EXF_INTERPOLATION))
180 c Shift wind stresses calculated at C-points to W/S points
181 fv(i,j,bi,bj) = exf_outscal_vstress*
182 & (vstress(i,j,bi,bj)+vstress(i,j-1,bi,bj))/2.*
183 & maskS(i,j,1,bi,bj)
184 #else
185 fv(i,j,bi,bj) = exf_outscal_vstress*vstress(i,j,bi,bj)
186 #endif
187 enddo
188 enddo
189
190 #ifdef SHORTWAVE_HEATING
191 c Short wave radiative flux.
192 do j = jmin,jmax
193 do i = imin,imax
194 qsw(i,j,bi,bj) = exf_outscal_swflux*swflux(i,j,bi,bj)
195 enddo
196 enddo
197 C net heat flux = heat flux (except SW) + SW flux.
198 do j = jmin,jmax
199 do i = imin,imax
200 qnet(i,j,bi,bj) = qnet(i,j,bi,bj) + qsw(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 c _EXCH_XY_R4( fu, mythid )
237 c _EXCH_XY_R4( fv, mythid )
238 CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid)
239 #ifdef SHORTWAVE_HEATING
240 _EXCH_XY_R4( qsw, mythid )
241 #endif
242 #ifdef ALLOW_CLIMSST_RELAXATION
243 _EXCH_XY_R4( sst, mythid )
244 #endif
245 #ifdef ALLOW_CLIMSSS_RELAXATION
246 _EXCH_XY_R4( sss, mythid )
247 #endif
248 #ifdef ATMOSPHERIC_LOADING
249 _EXCH_XY_R4( pload, mythid )
250 #endif
251
252 end

  ViewVC Help
Powered by ViewVC 1.1.22