/[MITgcm]/MITgcm/pkg/exf/exf_mapfields.F
ViewVC logotype

Annotation of /MITgcm/pkg/exf/exf_mapfields.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.24 - (hide annotations) (download)
Thu Dec 22 19:03:41 2011 UTC (12 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint63l, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k
Changes since 1.23: +8 -5 lines
remove/avoid un-used variables

1 jmc 1.24 C $Header: /u/gcmpack/MITgcm/pkg/exf/exf_mapfields.F,v 1.23 2009/04/28 18:15:33 jmc Exp $
2 jmc 1.17 C $Name: $
3 heimbach 1.1
4 edhill 1.7 #include "EXF_OPTIONS.h"
5 heimbach 1.1
6 heimbach 1.16 subroutine exf_mapfields( mytime, myiter, mythid )
7 heimbach 1.1
8     c ==================================================================
9 heimbach 1.2 c SUBROUTINE exf_mapfields
10 heimbach 1.1 c ==================================================================
11     c
12 dimitri 1.5 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 jmc 1.17 c EXF_FIELDS.h and FFIELDS.h for definitions of the various input
20 dimitri 1.5 c and output fields and for default unit and sign convetions.
21 heimbach 1.1 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 dimitri 1.5 c scaling factors exf_outscal_ust, exf_outscal_vst
38 dimitri 1.4 c
39 dimitri 1.5 c mods for pkg/seaice: menemenlis@jpl.nasa.gov 20-Dec-2002
40 dimitri 1.4 c
41 heimbach 1.1 c ==================================================================
42 heimbach 1.2 c SUBROUTINE exf_mapfields
43 heimbach 1.1 c ==================================================================
44    
45     implicit none
46    
47     c == global variables ==
48    
49     #include "EEPARAMS.h"
50     #include "SIZE.h"
51 heimbach 1.16 #include "PARAMS.h"
52 heimbach 1.1 #include "FFIELDS.h"
53 heimbach 1.8 #include "GRID.h"
54    
55 jmc 1.17 #include "EXF_PARAM.h"
56     #include "EXF_CONSTANTS.h"
57     #include "EXF_FIELDS.h"
58 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
59     # include "tamc.h"
60     # include "tamc_keys.h"
61     #endif
62 heimbach 1.1 c == routine arguments ==
63    
64     c mythid - thread number for this instance of the routine.
65    
66 jmc 1.20 _RL mytime
67     integer myiter
68 heimbach 1.1 integer mythid
69    
70     c == local variables ==
71    
72 jmc 1.24 INTEGER bi,bj
73     INTEGER i,j,ks
74 jmc 1.20 INTEGER imin, imax
75     INTEGER jmin, jmax
76     PARAMETER ( imin = 1-OLx , imax = sNx+OLx )
77     PARAMETER ( jmin = 1-OLy , jmax = sNy+OLy )
78 heimbach 1.1
79     c == end of interface ==
80    
81 jmc 1.24 C-- set surface level index:
82     ks = 1
83    
84 jmc 1.20 DO bj = myByLo(myThid),myByHi(myThid)
85     DO bi = myBxLo(myThid),myBxHi(myThid)
86 heimbach 1.2
87     #ifdef ALLOW_AUTODIFF_TAMC
88     act1 = bi - myBxLo(myThid)
89     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
90     act2 = bj - myByLo(myThid)
91     max2 = myByHi(myThid) - myByLo(myThid) + 1
92     act3 = myThid - 1
93     max3 = nTx*nTy
94     act4 = ikey_dynamics - 1
95     ikey = (act1 + 1) + act2*max1
96     & + act3*max1*max2
97     & + act4*max1*max2*max3
98     #endif /* ALLOW_AUTODIFF_TAMC */
99    
100 mlosch 1.21 c Heat flux.
101 heimbach 1.1 do j = jmin,jmax
102     do i = imin,imax
103 heimbach 1.6 qnet(i,j,bi,bj) = exf_outscal_hflux*hflux(i,j,bi,bj)
104 mlosch 1.21 enddo
105     enddo
106     if ( hfluxfile .EQ. ' ' ) then
107     do j = jmin,jmax
108     do i = imin,imax
109     qnet(i,j,bi,bj) = qnet(i,j,bi,bj) -
110 heimbach 1.16 & exf_outscal_hflux * ( hflux_exfremo_intercept +
111     & hflux_exfremo_slope*(mytime-starttime) )
112 heimbach 1.2 enddo
113 mlosch 1.21 enddo
114     endif
115 heimbach 1.2
116 mlosch 1.21 c Salt flux.
117 heimbach 1.2 do j = jmin,jmax
118     do i = imin,imax
119 jmc 1.22 EmPmR(i,j,bi,bj)= exf_outscal_sflux*sflux(i,j,bi,bj)
120     & *rhoConstFresh
121 mlosch 1.21 enddo
122     enddo
123     if ( sfluxfile .EQ. ' ' ) then
124     do j = jmin,jmax
125     do i = imin,imax
126 jmc 1.22 EmPmR(i,j,bi,bj) = EmPmR(i,j,bi,bj) - rhoConstFresh*
127 heimbach 1.16 & exf_outscal_sflux * ( sflux_exfremo_intercept +
128     & sflux_exfremo_slope*(mytime-starttime) )
129 heimbach 1.2 enddo
130 mlosch 1.21 enddo
131     endif
132 heimbach 1.1
133 heimbach 1.2 #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 heimbach 1.1 c Zonal wind stress.
139 mlosch 1.10 if (ustress(i,j,bi,bj).gt.windstressmax) then
140     ustress(i,j,bi,bj)=windstressmax
141 heimbach 1.2 endif
142     enddo
143     enddo
144     #ifdef ALLOW_AUTODIFF_TAMC
145     CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
146     #endif
147     do j = jmin,jmax
148     do i = imin,imax
149 mlosch 1.10 if (ustress(i,j,bi,bj).lt.-windstressmax) then
150     ustress(i,j,bi,bj)=-windstressmax
151 heimbach 1.2 endif
152     enddo
153     enddo
154 jmc 1.20 IF ( stressIsOnCgrid ) THEN
155     do j = jmin,jmax
156     do i = imin+1,imax
157     fu(i,j,bi,bj) = exf_outscal_ustress*ustress(i,j,bi,bj)
158     enddo
159     enddo
160     ELSE
161     do j = jmin,jmax
162 heimbach 1.8 do i = imin+1,imax
163 jmc 1.20 c Shift wind stresses calculated at Grid-center to W/S points
164 heimbach 1.8 fu(i,j,bi,bj) = exf_outscal_ustress*
165 jmc 1.20 & (ustress(i,j,bi,bj)+ustress(i-1,j,bi,bj))
166 jmc 1.24 & *exf_half*maskW(i,j,ks,bi,bj)
167 heimbach 1.2 enddo
168 jmc 1.20 enddo
169     ENDIF
170 heimbach 1.1
171 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
172     CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
173     #endif
174     do j = jmin,jmax
175     do i = imin,imax
176 heimbach 1.1 c Meridional wind stress.
177 mlosch 1.10 if (vstress(i,j,bi,bj).gt.windstressmax) then
178     vstress(i,j,bi,bj)=windstressmax
179 heimbach 1.2 endif
180     enddo
181     enddo
182     #ifdef ALLOW_AUTODIFF_TAMC
183     CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
184     #endif
185     do j = jmin,jmax
186     do i = imin,imax
187 mlosch 1.10 if (vstress(i,j,bi,bj).lt.-windstressmax) then
188     vstress(i,j,bi,bj)=-windstressmax
189 heimbach 1.2 endif
190     enddo
191     enddo
192 jmc 1.20 IF ( stressIsOnCgrid ) THEN
193     do j = jmin+1,jmax
194     do i = imin,imax
195     fv(i,j,bi,bj) = exf_outscal_vstress*vstress(i,j,bi,bj)
196     enddo
197     enddo
198     ELSE
199     do j = jmin+1,jmax
200 heimbach 1.2 do i = imin,imax
201 heimbach 1.8 c Shift wind stresses calculated at C-points to W/S points
202     fv(i,j,bi,bj) = exf_outscal_vstress*
203 jmc 1.20 & (vstress(i,j,bi,bj)+vstress(i,j-1,bi,bj))
204 jmc 1.24 & *exf_half*maskS(i,j,ks,bi,bj)
205 heimbach 1.2 enddo
206 jmc 1.20 enddo
207     ENDIF
208 heimbach 1.1
209 dimitri 1.5 #ifdef SHORTWAVE_HEATING
210 heimbach 1.1 c Short wave radiative flux.
211 heimbach 1.2 do j = jmin,jmax
212     do i = imin,imax
213 heimbach 1.6 qsw(i,j,bi,bj) = exf_outscal_swflux*swflux(i,j,bi,bj)
214 heimbach 1.2 enddo
215     enddo
216 heimbach 1.1 #endif
217    
218     #ifdef ALLOW_CLIMSST_RELAXATION
219 heimbach 1.2 do j = jmin,jmax
220     do i = imin,imax
221 dimitri 1.5 sst(i,j,bi,bj) = exf_outscal_sst*climsst(i,j,bi,bj)
222 heimbach 1.2 enddo
223     enddo
224 heimbach 1.1 #endif
225    
226     #ifdef ALLOW_CLIMSSS_RELAXATION
227 heimbach 1.2 do j = jmin,jmax
228     do i = imin,imax
229 dimitri 1.5 sss(i,j,bi,bj) = exf_outscal_sss*climsss(i,j,bi,bj)
230 heimbach 1.2 enddo
231     enddo
232 heimbach 1.1 #endif
233    
234 heimbach 1.2 #ifdef ATMOSPHERIC_LOADING
235     do j = jmin,jmax
236     do i = imin,imax
237 dimitri 1.5 pload(i,j,bi,bj)=exf_outscal_apressure*apressure(i,j,bi,bj)
238 heimbach 1.1 enddo
239     enddo
240 heimbach 1.2 #endif
241    
242 jmc 1.20 ENDDO
243     ENDDO
244 heimbach 1.1
245     c Update the tile edges.
246    
247 jmc 1.23 _EXCH_XY_RS( qnet, mythid )
248     _EXCH_XY_RS( empmr, mythid )
249 cheisey 1.3 CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid)
250 dimitri 1.5 #ifdef SHORTWAVE_HEATING
251 jmc 1.23 _EXCH_XY_RS( qsw, mythid )
252 heimbach 1.1 #endif
253     #ifdef ALLOW_CLIMSST_RELAXATION
254 jmc 1.23 _EXCH_XY_RS( sst, mythid )
255 heimbach 1.1 #endif
256     #ifdef ALLOW_CLIMSSS_RELAXATION
257 jmc 1.23 _EXCH_XY_RS( sss, mythid )
258 heimbach 1.2 #endif
259     #ifdef ATMOSPHERIC_LOADING
260 jmc 1.23 _EXCH_XY_RS( pload, mythid )
261 heimbach 1.1 #endif
262    
263 jmc 1.20 RETURN
264     END

  ViewVC Help
Powered by ViewVC 1.1.22