/[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.3 - (show annotations) (download)
Tue Dec 17 19:47:41 2002 UTC (21 years, 6 months ago) by cheisey
Branch: MAIN
CVS Tags: checkpoint47e_post
Branch point for: branch-exfmods-curt
Changes since 1.2: +4 -3 lines
Updating the exchange routines for fu, fv, ustress, vstress
for the cubed-sphere model.

1 c $Header: /u/u0/gcmpack/MITgcm/pkg/exf/exf_mapfields.F,v 1.2 2002/11/12 20:34:41 heimbach Exp $
2
3 #include "EXF_CPPOPTIONS.h"
4
5
6 subroutine exf_mapfields( mythid )
7
8 c ==================================================================
9 c SUBROUTINE exf_mapfields
10 c ==================================================================
11 c
12 c o Map the external forcing fields on the ocean model arrays. This
13 c routine is included to separate the ocean state estimation tool
14 c as much as possible from the ocean model. Unit conversion factors
15 c are to be set by the user.
16 c
17 c The units have to be such that the individual forcing record has
18 c units equal to [quantity/s]. See the header file *FFIELDS.h* of
19 c the MITgcmuv.
20 c
21 c Required units such that no scaling has to be applied:
22 c
23 c heat flux: input file W/m^2
24 c salt flux: input file m/s
25 c zonal wind stress: input file N/m^2
26 c merid. wind stress: input file N/m^2
27 c
28 c To allow for such unit conversions this routine contains scaling
29 c factors scal_quantity.
30 c
31 c started: Christian Eckert eckert@mit.edu 09-Aug-1999
32 c
33 c changed: Christian Eckert eckert@mit.edu 11-Jan-2000
34 c
35 c - Restructured the code in order to create a package
36 c for the MITgcmUV.
37 c
38 c Christian Eckert eckert@mit.edu 12-Feb-2000
39 c
40 c - Changed Routine names (package prefix: exf_)
41 c
42 c Patrick Heimbach, heimbach@mit.edu 06-May-2000
43 c
44 c - added and changed CPP flag structure for
45 c ALLOW_BULKFORMULAE, ALLOW_ATM_TEMP
46 c
47 c Patrick Heimbach, heimbach@mit.edu 23-May-2000
48 c
49 c - sign change of ustress/vstress incorporated into
50 c scaling factors scal_ust, scal_vst
51 c
52 c ==================================================================
53 c SUBROUTINE exf_mapfields
54 c ==================================================================
55
56 implicit none
57
58 c == global variables ==
59
60 #include "EEPARAMS.h"
61 #include "SIZE.h"
62 #include "FFIELDS.h"
63 #include "exf_param.h"
64 #include "exf_constants.h"
65 #include "exf_fields.h"
66 #include "exf_clim_fields.h"
67 #ifdef ALLOW_AUTODIFF_TAMC
68 # include "tamc.h"
69 # include "tamc_keys.h"
70 #endif
71 c == routine arguments ==
72
73 c mythid - thread number for this instance of the routine.
74
75 integer mythid
76
77 c == local variables ==
78
79 integer bi,bj
80 integer i,j
81 integer jtlo
82 integer jthi
83 integer itlo
84 integer ithi
85 integer jmin
86 integer jmax
87 integer imin
88 integer imax
89
90 c == end of interface ==
91
92 jtlo = mybylo(mythid)
93 jthi = mybyhi(mythid)
94 itlo = mybxlo(mythid)
95 ithi = mybxhi(mythid)
96 jmin = 1-oly
97 jmax = sny+oly
98 imin = 1-olx
99 imax = snx+olx
100
101 do bj = jtlo,jthi
102 do bi = itlo,ithi
103
104 #ifdef ALLOW_AUTODIFF_TAMC
105 act1 = bi - myBxLo(myThid)
106 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
107 act2 = bj - myByLo(myThid)
108 max2 = myByHi(myThid) - myByLo(myThid) + 1
109 act3 = myThid - 1
110 max3 = nTx*nTy
111 act4 = ikey_dynamics - 1
112 ikey = (act1 + 1) + act2*max1
113 & + act3*max1*max2
114 & + act4*max1*max2*max3
115 #endif /* ALLOW_AUTODIFF_TAMC */
116
117 do j = jmin,jmax
118 do i = imin,imax
119 c Heat flux.
120 qnet(i,j,bi,bj) = scal_hfl*hflux(i,j,bi,bj)
121 enddo
122 enddo
123
124
125 do j = jmin,jmax
126 do i = imin,imax
127 c Salt flux.
128 #if (defined (ALLOW_BULKFORMULAE) && defined (ALLOW_ATM_TEMP))
129 empmr(i,j,bi,bj) = scal_prc*precip(i,j,bi,bj)
130 #else
131 empmr(i,j,bi,bj) = scal_sfl*sflux(i,j,bi,bj)
132 #endif
133 enddo
134 enddo
135
136 #ifdef ALLOW_AUTODIFF_TAMC
137 CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
138 #endif
139 do j = jmin,jmax
140 do i = imin,imax
141 c Zonal wind stress.
142 if (ustress(i,j,bi,bj).gt.2.0D0) then
143 ustress(i,j,bi,bj)=2.0D0
144 endif
145 enddo
146 enddo
147 #ifdef ALLOW_AUTODIFF_TAMC
148 CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
149 #endif
150 do j = jmin,jmax
151 do i = imin,imax
152 if (ustress(i,j,bi,bj).lt.-2.0D0) then
153 ustress(i,j,bi,bj)=-2.0D0
154 endif
155 enddo
156 enddo
157 do j = jmin,jmax
158 do i = imin,imax
159 fu(i,j,bi,bj) = scal_ust*ustress(i,j,bi,bj)
160 enddo
161 enddo
162
163 #ifdef ALLOW_AUTODIFF_TAMC
164 CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
165 #endif
166 do j = jmin,jmax
167 do i = imin,imax
168 c Meridional wind stress.
169 if (vstress(i,j,bi,bj).gt.2.0D0) then
170 vstress(i,j,bi,bj)=2.0D0
171 endif
172 enddo
173 enddo
174 #ifdef ALLOW_AUTODIFF_TAMC
175 CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
176 #endif
177 do j = jmin,jmax
178 do i = imin,imax
179 if (vstress(i,j,bi,bj).lt.-2.0D0) then
180 vstress(i,j,bi,bj)=-2.0D0
181 endif
182 enddo
183 enddo
184 do j = jmin,jmax
185 do i = imin,imax
186 fv(i,j,bi,bj) = scal_vst*vstress(i,j,bi,bj)
187 enddo
188 enddo
189
190 #ifdef ALLOW_KPP || \
191 (defined (ALLOW_BULKFORMULAE) && defined (ALLOW_ATM_TEMP)))
192 c Short wave radiative flux.
193 do j = jmin,jmax
194 do i = imin,imax
195 qsw(i,j,bi,bj) = scal_swf*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) = scal_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) = scal_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) = scal_apressure*apressure(i,j,bi,bj)
220 enddo
221 enddo
222 #endif
223
224
225 enddo
226 enddo
227
228 c Update the tile edges.
229
230 _EXCH_XY_R4( qnet, mythid )
231 _EXCH_XY_R4( empmr, mythid )
232 c _EXCH_XY_R4( fu, mythid )
233 c _EXCH_XY_R4( fv, mythid )
234 CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid)
235 #ifdef ALLOW_KPP
236 _EXCH_XY_R4( qsw, mythid )
237 #endif
238 #ifdef ALLOW_CLIMSST_RELAXATION
239 _EXCH_XY_R4( sst, mythid )
240 #endif
241 #ifdef ALLOW_CLIMSSS_RELAXATION
242 _EXCH_XY_R4( sss, mythid )
243 #endif
244 #ifdef ATMOSPHERIC_LOADING
245 _EXCH_XY_R4( pload, mythid )
246 #endif
247
248 end

  ViewVC Help
Powered by ViewVC 1.1.22