/[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.16 - (show annotations) (download)
Thu Mar 2 15:30:11 2006 UTC (18 years, 2 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58b_post, checkpoint58f_post, checkpoint58d_post, checkpoint58y_post, checkpoint58t_post, checkpoint58m_post, checkpoint58w_post, checkpoint58o_post, checkpoint58p_post, checkpoint58q_post, checkpoint58e_post, checkpoint58r_post, checkpoint58n_post, checkpoint58k_post, checkpoint58v_post, checkpoint58l_post, checkpoint58g_post, checkpoint58x_post, checkpoint58h_post, checkpoint58j_post, checkpoint58i_post, checkpoint58c_post, checkpoint58u_post, checkpoint58s_post
Changes since 1.15: +13 -2 lines
o Clean exf namelist
o Update trend removal code

1 c $Header: /u/gcmpack/MITgcm/pkg/exf/exf_mapfields.F,v 1.15 2006/01/13 16:48:02 heimbach Exp $
2
3 #include "EXF_OPTIONS.h"
4
5 subroutine exf_mapfields( mytime, myiter, 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 "PARAMS.h"
51 #include "FFIELDS.h"
52 #include "GRID.h"
53
54 #include "exf_param.h"
55 #include "exf_constants.h"
56 #include "exf_fields.h"
57 #include "exf_clim_param.h"
58 #include "exf_clim_fields.h"
59 #ifdef ALLOW_AUTODIFF_TAMC
60 # include "tamc.h"
61 # include "tamc_keys.h"
62 #endif
63 c == routine arguments ==
64
65 c mythid - thread number for this instance of the routine.
66
67 integer mythid
68 integer myiter
69 _RL mytime
70
71 c == local variables ==
72
73 integer bi,bj
74 integer i,j,k
75 integer jtlo
76 integer jthi
77 integer itlo
78 integer ithi
79 integer jmin
80 integer jmax
81 integer imin
82 integer imax
83
84 c == end of interface ==
85
86 jtlo = mybylo(mythid)
87 jthi = mybyhi(mythid)
88 itlo = mybxlo(mythid)
89 ithi = mybxhi(mythid)
90 jmin = 1-oly
91 jmax = sny+oly
92 imin = 1-olx
93 imax = snx+olx
94
95 do bj = jtlo,jthi
96 do bi = itlo,ithi
97
98 #ifdef ALLOW_AUTODIFF_TAMC
99 act1 = bi - myBxLo(myThid)
100 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
101 act2 = bj - myByLo(myThid)
102 max2 = myByHi(myThid) - myByLo(myThid) + 1
103 act3 = myThid - 1
104 max3 = nTx*nTy
105 act4 = ikey_dynamics - 1
106 ikey = (act1 + 1) + act2*max1
107 & + act3*max1*max2
108 & + act4*max1*max2*max3
109 #endif /* ALLOW_AUTODIFF_TAMC */
110
111 do j = jmin,jmax
112 do i = imin,imax
113 c Heat flux.
114 qnet(i,j,bi,bj) = exf_outscal_hflux*hflux(i,j,bi,bj)
115 if ( hfluxfile .EQ. ' ' )
116 & qnet(i,j,bi,bj) = qnet(i,j,bi,bj) -
117 & exf_outscal_hflux * ( hflux_exfremo_intercept +
118 & hflux_exfremo_slope*(mytime-starttime) )
119 enddo
120 enddo
121
122
123 do j = jmin,jmax
124 do i = imin,imax
125 c Salt flux.
126 empmr(i,j,bi,bj)= exf_outscal_sflux*sflux(i,j,bi,bj)
127 if ( sfluxfile .EQ. ' ' )
128 & empmr(i,j,bi,bj) = empmr(i,j,bi,bj) -
129 & exf_outscal_sflux * ( sflux_exfremo_intercept +
130 & sflux_exfremo_slope*(mytime-starttime) )
131 enddo
132 enddo
133
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 c Zonal wind stress.
140 if (ustress(i,j,bi,bj).gt.windstressmax) then
141 ustress(i,j,bi,bj)=windstressmax
142 endif
143 enddo
144 enddo
145 #ifdef ALLOW_AUTODIFF_TAMC
146 CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
147 #endif
148 do j = jmin,jmax
149 do i = imin,imax
150 if (ustress(i,j,bi,bj).lt.-windstressmax) then
151 ustress(i,j,bi,bj)=-windstressmax
152 endif
153 enddo
154 enddo
155 do j = jmin,jmax
156 do i = imin+1,imax
157 #if (defined (ALLOW_BULKFORMULAE) || defined (USE_EXF_INTERPOLATION))
158 c Shift wind stresses calculated at C-points to W/S points
159 fu(i,j,bi,bj) = exf_outscal_ustress*
160 & (ustress(i,j,bi,bj)+ustress(i-1,j,bi,bj))/2.*
161 & maskW(i,j,1,bi,bj)
162 #else
163 fu(i,j,bi,bj) = exf_outscal_ustress*ustress(i,j,bi,bj)
164 #endif
165 enddo
166 enddo
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 do j = jmin+1,jmax
190 do i = imin,imax
191 #if (defined (ALLOW_BULKFORMULAE) || defined (USE_EXF_INTERPOLATION))
192 c Shift wind stresses calculated at C-points to W/S points
193 fv(i,j,bi,bj) = exf_outscal_vstress*
194 & (vstress(i,j,bi,bj)+vstress(i,j-1,bi,bj))/2.*
195 & maskS(i,j,1,bi,bj)
196 #else
197 fv(i,j,bi,bj) = exf_outscal_vstress*vstress(i,j,bi,bj)
198 #endif
199 enddo
200 enddo
201
202 #ifdef SHORTWAVE_HEATING
203 c Short wave radiative flux.
204 do j = jmin,jmax
205 do i = imin,imax
206 qsw(i,j,bi,bj) = exf_outscal_swflux*swflux(i,j,bi,bj)
207 enddo
208 enddo
209 #endif
210
211 #ifdef ALLOW_CLIMSST_RELAXATION
212 do j = jmin,jmax
213 do i = imin,imax
214 sst(i,j,bi,bj) = exf_outscal_sst*climsst(i,j,bi,bj)
215 enddo
216 enddo
217 #endif
218
219 #ifdef ALLOW_CLIMSSS_RELAXATION
220 do j = jmin,jmax
221 do i = imin,imax
222 sss(i,j,bi,bj) = exf_outscal_sss*climsss(i,j,bi,bj)
223 enddo
224 enddo
225 #endif
226
227 #ifdef ATMOSPHERIC_LOADING
228 do j = jmin,jmax
229 do i = imin,imax
230 pload(i,j,bi,bj)=exf_outscal_apressure*apressure(i,j,bi,bj)
231 enddo
232 enddo
233 #endif
234
235 enddo
236 enddo
237
238 c Update the tile edges.
239
240 _EXCH_XY_R4( qnet, mythid )
241 _EXCH_XY_R4( empmr, mythid )
242 c _EXCH_XY_R4( fu, mythid )
243 c _EXCH_XY_R4( fv, mythid )
244 CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid)
245 #ifdef SHORTWAVE_HEATING
246 _EXCH_XY_R4( qsw, mythid )
247 #endif
248 #ifdef ALLOW_CLIMSST_RELAXATION
249 _EXCH_XY_R4( sst, mythid )
250 #endif
251 #ifdef ALLOW_CLIMSSS_RELAXATION
252 _EXCH_XY_R4( sss, mythid )
253 #endif
254 #ifdef ATMOSPHERIC_LOADING
255 _EXCH_XY_R4( pload, mythid )
256 #endif
257
258 end

  ViewVC Help
Powered by ViewVC 1.1.22