/[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.21 - (hide annotations) (download)
Thu Sep 27 09:41:13 2007 UTC (16 years, 8 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint59h
Changes since 1.20: +19 -10 lines
improve vectorization by moving if-statement out of loops

1 mlosch 1.21 C $Header: /u/gcmpack/MITgcm/pkg/exf/exf_mapfields.F,v 1.20 2007/05/14 19:34:57 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     integer bi,bj
73 heimbach 1.13 integer i,j,k
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.20 DO bj = myByLo(myThid),myByHi(myThid)
82     DO bi = myBxLo(myThid),myBxHi(myThid)
83 heimbach 1.2
84     #ifdef ALLOW_AUTODIFF_TAMC
85     act1 = bi - myBxLo(myThid)
86     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
87     act2 = bj - myByLo(myThid)
88     max2 = myByHi(myThid) - myByLo(myThid) + 1
89     act3 = myThid - 1
90     max3 = nTx*nTy
91     act4 = ikey_dynamics - 1
92     ikey = (act1 + 1) + act2*max1
93     & + act3*max1*max2
94     & + act4*max1*max2*max3
95     #endif /* ALLOW_AUTODIFF_TAMC */
96    
97 mlosch 1.21 c Heat flux.
98 heimbach 1.1 do j = jmin,jmax
99     do i = imin,imax
100 heimbach 1.6 qnet(i,j,bi,bj) = exf_outscal_hflux*hflux(i,j,bi,bj)
101 mlosch 1.21 enddo
102     enddo
103     if ( hfluxfile .EQ. ' ' ) then
104     do j = jmin,jmax
105     do i = imin,imax
106     qnet(i,j,bi,bj) = qnet(i,j,bi,bj) -
107 heimbach 1.16 & exf_outscal_hflux * ( hflux_exfremo_intercept +
108     & hflux_exfremo_slope*(mytime-starttime) )
109 heimbach 1.2 enddo
110 mlosch 1.21 enddo
111     endif
112 heimbach 1.2
113 mlosch 1.21 c Salt flux.
114 heimbach 1.2 do j = jmin,jmax
115     do i = imin,imax
116 heimbach 1.6 empmr(i,j,bi,bj)= exf_outscal_sflux*sflux(i,j,bi,bj)
117 mlosch 1.21 enddo
118     enddo
119     if ( sfluxfile .EQ. ' ' ) then
120     do j = jmin,jmax
121     do i = imin,imax
122     empmr(i,j,bi,bj) = empmr(i,j,bi,bj) -
123 heimbach 1.16 & exf_outscal_sflux * ( sflux_exfremo_intercept +
124     & sflux_exfremo_slope*(mytime-starttime) )
125 heimbach 1.2 enddo
126 mlosch 1.21 enddo
127     endif
128 heimbach 1.1
129 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
130     CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
131     #endif
132     do j = jmin,jmax
133     do i = imin,imax
134 heimbach 1.1 c Zonal wind stress.
135 mlosch 1.10 if (ustress(i,j,bi,bj).gt.windstressmax) then
136     ustress(i,j,bi,bj)=windstressmax
137 heimbach 1.2 endif
138     enddo
139     enddo
140     #ifdef ALLOW_AUTODIFF_TAMC
141     CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
142     #endif
143     do j = jmin,jmax
144     do i = imin,imax
145 mlosch 1.10 if (ustress(i,j,bi,bj).lt.-windstressmax) then
146     ustress(i,j,bi,bj)=-windstressmax
147 heimbach 1.2 endif
148     enddo
149     enddo
150 jmc 1.20 IF ( stressIsOnCgrid ) THEN
151     do j = jmin,jmax
152     do i = imin+1,imax
153     fu(i,j,bi,bj) = exf_outscal_ustress*ustress(i,j,bi,bj)
154     enddo
155     enddo
156     ELSE
157     do j = jmin,jmax
158 heimbach 1.8 do i = imin+1,imax
159 jmc 1.20 c Shift wind stresses calculated at Grid-center to W/S points
160 heimbach 1.8 fu(i,j,bi,bj) = exf_outscal_ustress*
161 jmc 1.20 & (ustress(i,j,bi,bj)+ustress(i-1,j,bi,bj))
162     & *exf_half*maskW(i,j,1,bi,bj)
163 heimbach 1.2 enddo
164 jmc 1.20 enddo
165     ENDIF
166 heimbach 1.1
167 heimbach 1.2 #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 heimbach 1.1 c Meridional wind stress.
173 mlosch 1.10 if (vstress(i,j,bi,bj).gt.windstressmax) then
174     vstress(i,j,bi,bj)=windstressmax
175 heimbach 1.2 endif
176     enddo
177     enddo
178     #ifdef ALLOW_AUTODIFF_TAMC
179     CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
180     #endif
181     do j = jmin,jmax
182     do i = imin,imax
183 mlosch 1.10 if (vstress(i,j,bi,bj).lt.-windstressmax) then
184     vstress(i,j,bi,bj)=-windstressmax
185 heimbach 1.2 endif
186     enddo
187     enddo
188 jmc 1.20 IF ( stressIsOnCgrid ) THEN
189     do j = jmin+1,jmax
190     do i = imin,imax
191     fv(i,j,bi,bj) = exf_outscal_vstress*vstress(i,j,bi,bj)
192     enddo
193     enddo
194     ELSE
195     do j = jmin+1,jmax
196 heimbach 1.2 do i = imin,imax
197 heimbach 1.8 c Shift wind stresses calculated at C-points to W/S points
198     fv(i,j,bi,bj) = exf_outscal_vstress*
199 jmc 1.20 & (vstress(i,j,bi,bj)+vstress(i,j-1,bi,bj))
200     & *exf_half*maskS(i,j,1,bi,bj)
201 heimbach 1.2 enddo
202 jmc 1.20 enddo
203     ENDIF
204 heimbach 1.1
205 dimitri 1.5 #ifdef SHORTWAVE_HEATING
206 heimbach 1.1 c Short wave radiative flux.
207 heimbach 1.2 do j = jmin,jmax
208     do i = imin,imax
209 heimbach 1.6 qsw(i,j,bi,bj) = exf_outscal_swflux*swflux(i,j,bi,bj)
210 heimbach 1.2 enddo
211     enddo
212 heimbach 1.1 #endif
213    
214     #ifdef ALLOW_CLIMSST_RELAXATION
215 heimbach 1.2 do j = jmin,jmax
216     do i = imin,imax
217 dimitri 1.5 sst(i,j,bi,bj) = exf_outscal_sst*climsst(i,j,bi,bj)
218 heimbach 1.2 enddo
219     enddo
220 heimbach 1.1 #endif
221    
222     #ifdef ALLOW_CLIMSSS_RELAXATION
223 heimbach 1.2 do j = jmin,jmax
224     do i = imin,imax
225 dimitri 1.5 sss(i,j,bi,bj) = exf_outscal_sss*climsss(i,j,bi,bj)
226 heimbach 1.2 enddo
227     enddo
228 heimbach 1.1 #endif
229    
230 heimbach 1.2 #ifdef ATMOSPHERIC_LOADING
231     do j = jmin,jmax
232     do i = imin,imax
233 dimitri 1.5 pload(i,j,bi,bj)=exf_outscal_apressure*apressure(i,j,bi,bj)
234 heimbach 1.1 enddo
235     enddo
236 heimbach 1.2 #endif
237    
238 jmc 1.20 ENDDO
239     ENDDO
240 heimbach 1.1
241     c Update the tile edges.
242    
243     _EXCH_XY_R4( qnet, mythid )
244     _EXCH_XY_R4( empmr, mythid )
245 cheisey 1.3 CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid)
246 dimitri 1.5 #ifdef SHORTWAVE_HEATING
247 heimbach 1.1 _EXCH_XY_R4( qsw, mythid )
248     #endif
249     #ifdef ALLOW_CLIMSST_RELAXATION
250     _EXCH_XY_R4( sst, mythid )
251     #endif
252     #ifdef ALLOW_CLIMSSS_RELAXATION
253     _EXCH_XY_R4( sss, mythid )
254 heimbach 1.2 #endif
255     #ifdef ATMOSPHERIC_LOADING
256     _EXCH_XY_R4( pload, mythid )
257 heimbach 1.1 #endif
258    
259 jmc 1.20 RETURN
260     END

  ViewVC Help
Powered by ViewVC 1.1.22