/[MITgcm]/MITgcm/pkg/ctrl/ctrl_map_ini.F
ViewVC logotype

Annotation of /MITgcm/pkg/ctrl/ctrl_map_ini.F

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


Revision 1.15 - (hide annotations) (download)
Mon Feb 28 17:29:38 2005 UTC (19 years, 3 months ago) by heimbach
Branch: MAIN
Changes since 1.14: +53 -2 lines
Adding eddy stress controls a la Ferreira et al.

1 heimbach 1.15 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_map_ini.F,v 1.14 2004/11/16 05:42:12 heimbach Exp $
2 heimbach 1.1
3     #include "CTRL_CPPOPTIONS.h"
4    
5 heimbach 1.5 CBOP
6     C !ROUTINE: ctrl_map_ini
7     C !INTERFACE:
8     subroutine ctrl_map_ini( mythid )
9    
10     C !DESCRIPTION: \bv
11     c *=================================================================
12     c | SUBROUTINE ctrl_map_ini
13     c | Add the temperature, salinity, and diffusivity parts of the
14     c | control vector to the model state and update the tile halos.
15     c | The control vector is defined in the header file "ctrl.h".
16     c *=================================================================
17     C \ev
18 heimbach 1.1
19 heimbach 1.5 C !USES:
20 heimbach 1.1 implicit none
21    
22     c == global variables ==
23 heimbach 1.6 #include "SIZE.h"
24 heimbach 1.1 #include "EEPARAMS.h"
25 heimbach 1.6 #include "PARAMS.h"
26 heimbach 1.15 #include "GRID.h"
27 heimbach 1.1 #include "DYNVARS.h"
28 heimbach 1.15 #include "FFIELDS.h"
29 heimbach 1.1 #include "ctrl.h"
30     #include "ctrl_dummy.h"
31 heimbach 1.2 #include "optim.h"
32 heimbach 1.13 #ifdef ALLOW_PTRACERS
33     # include "PTRACERS_SIZE.h"
34     # include "PTRACERS.h"
35     #endif
36 heimbach 1.12 #ifdef ALLOW_ECCO
37     # include "ecco_cost.h"
38     #endif
39 heimbach 1.1
40 heimbach 1.5 C !INPUT/OUTPUT PARAMETERS:
41 heimbach 1.1 c == routine arguments ==
42     integer mythid
43    
44 heimbach 1.5 C !LOCAL VARIABLES:
45 heimbach 1.1 c == local variables ==
46    
47     integer bi,bj
48     integer i,j,k
49     integer itlo,ithi
50     integer jtlo,jthi
51     integer jmin,jmax
52     integer imin,imax
53     integer il
54    
55     logical equal
56     logical doglobalread
57     logical ladinit
58    
59     character*( 80) fnametheta
60     character*( 80) fnamesalt
61 heimbach 1.2 character*( 80) fnametr1
62 heimbach 1.3 character*( 80) fnamediffkr
63     character*( 80) fnamekapgm
64 heimbach 1.6 character*( 80) fnameefluxy
65     character*( 80) fnameefluxp
66 heimbach 1.7 character*( 80) fnamebottomdrag
67 heimbach 1.14 character*( 80) fnamesss
68     character*( 80) fnamesst
69 heimbach 1.15 character*( 80) fnameedtaux
70     character*( 80) fnameedtauy
71 heimbach 1.1
72 heimbach 1.5 _RL fac
73 heimbach 1.12 _RL tmptest
74 heimbach 1.5
75 heimbach 1.1 c == external ==
76     integer ilnblnk
77     external ilnblnk
78    
79     c == end of interface ==
80 heimbach 1.5 CEOP
81 heimbach 1.1
82     jtlo = mybylo(mythid)
83     jthi = mybyhi(mythid)
84     itlo = mybxlo(mythid)
85     ithi = mybxhi(mythid)
86 heimbach 1.8 jmin = 1
87     jmax = sny
88     imin = 1
89     imax = snx
90 heimbach 1.1
91     doglobalread = .false.
92     ladinit = .false.
93    
94     equal = .true.
95    
96     if ( equal ) then
97     fac = 1. _d 0
98     else
99     fac = 0. _d 0
100     endif
101    
102     #ifdef ALLOW_THETA0_CONTROL
103     c-- Temperature field.
104     il=ilnblnk( xx_theta_file )
105     write(fnametheta(1:80),'(2a,i10.10)')
106     & xx_theta_file(1:il),'.',optimcycle
107 heimbach 1.11 call active_read_xyz_loc( fnametheta, tmpfld3d, 1,
108 heimbach 1.1 & doglobalread, ladinit, optimcycle,
109     & mythid, xx_theta_dummy )
110    
111     do bj = jtlo,jthi
112     do bi = itlo,ithi
113     do k = 1,nr
114     do j = jmin,jmax
115     do i = imin,imax
116 heimbach 1.12 #ifdef ALLOW_ECCO
117     IF (abs(tmpfld3d(i,j,k,bi,bj)).gt.
118     $ 2.0/sqrt(wtheta(k,bi,bj)))
119     $ tmpfld3d(i,j,k,bi,bj)=
120     $ sign(2.0/sqrt(wtheta(k,bi,bj)),tmpfld3d(i,j,k,bi,bj))
121     #endif
122 heimbach 1.1 theta(i,j,k,bi,bj) = theta(i,j,k,bi,bj) +
123     & fac*tmpfld3d(i,j,k,bi,bj)
124 heimbach 1.9 if(theta(i,j,k,bi,bj).lt.-2.0)
125 heimbach 1.12 & theta(i,j,k,bi,bj)= -2.0
126 heimbach 1.1 enddo
127     enddo
128     enddo
129     enddo
130     enddo
131 heimbach 1.12
132 heimbach 1.1 #endif
133    
134     #ifdef ALLOW_SALT0_CONTROL
135     c-- Temperature field.
136     il=ilnblnk( xx_salt_file )
137     write(fnamesalt(1:80),'(2a,i10.10)')
138     & xx_salt_file(1:il),'.',optimcycle
139 heimbach 1.11 call active_read_xyz_loc( fnamesalt, tmpfld3d, 1,
140 heimbach 1.1 & doglobalread, ladinit, optimcycle,
141     & mythid, xx_salt_dummy )
142    
143     do bj = jtlo,jthi
144     do bi = itlo,ithi
145     do k = 1,nr
146     do j = jmin,jmax
147     do i = imin,imax
148 heimbach 1.12 #ifdef ALLOW_ECCO
149     IF (abs(tmpfld3d(i,j,k,bi,bj)).gt.
150     $ 2.0/sqrt(wsalt(k,bi,bj)))
151     $ tmpfld3d(i,j,k,bi,bj)=
152     $ sign(2.0/sqrt(wsalt(k,bi,bj)),tmpfld3d(i,j,k,bi,bj))
153     #endif
154 heimbach 1.1 salt(i,j,k,bi,bj) = salt(i,j,k,bi,bj) +
155     & fac*tmpfld3d(i,j,k,bi,bj)
156 heimbach 1.12
157 heimbach 1.1 enddo
158     enddo
159     enddo
160     enddo
161 heimbach 1.12 enddO
162 heimbach 1.1 #endif
163    
164 heimbach 1.2 #ifdef ALLOW_TR10_CONTROL
165 heimbach 1.13 #ifdef ALLOW_PTRACERS
166 heimbach 1.2 c-- Temperature field.
167     il=ilnblnk( xx_tr1_file )
168     write(fnametr1(1:80),'(2a,i10.10)')
169     & xx_tr1_file(1:il),'.',optimcycle
170 heimbach 1.11 call active_read_xyz_loc( fnametr1, tmpfld3d, 1,
171 heimbach 1.2 & doglobalread, ladinit, optimcycle,
172     & mythid, xx_tr1_dummy )
173    
174     do bj = jtlo,jthi
175     do bi = itlo,ithi
176     do k = 1,nr
177     do j = jmin,jmax
178     do i = imin,imax
179 heimbach 1.13 ptracer(i,j,k,bi,bj,1) = ptracer(i,j,k,bi,bj,1) +
180 heimbach 1.2 & fac*tmpfld3d(i,j,k,bi,bj)
181     enddo
182     enddo
183     enddo
184     enddo
185     enddo
186     #endif
187 heimbach 1.13 #endif
188 heimbach 1.2
189 heimbach 1.14 #ifdef ALLOW_SST0_CONTROL
190     c-- sst0.
191     il=ilnblnk( xx_sst_file )
192     write(fnamesst(1:80),'(2a,i10.10)')
193     & xx_sst_file(1:il),'.',optimcycle
194     call active_read_xy_loc ( fnamesst, tmpfld2d, 1,
195     & doglobalread, ladinit, optimcycle,
196     & mythid, xx_sst_dummy )
197     do bj = jtlo,jthi
198     do bi = itlo,ithi
199     do j = jmin,jmax
200     do i = imin,imax
201     cph sst(i,j,bi,bj) = sst(i,j,bi,bj) + tmpfld2d(i,j,bi,bj)
202     theta(i,j,1,bi,bj) = theta(i,j,1,bi,bj)
203     & + tmpfld2d(i,j,bi,bj)
204     enddo
205     enddo
206     enddo
207     enddo
208     #endif
209    
210     #ifdef ALLOW_SSS0_CONTROL
211     c-- sss0.
212     il=ilnblnk( xx_sss_file )
213     write(fnamesss(1:80),'(2a,i10.10)')
214     & xx_sss_file(1:il),'.',optimcycle
215     call active_read_xy_loc ( fnamesss, tmpfld2d, 1,
216     & doglobalread, ladinit, optimcycle,
217     & mythid, xx_sss_dummy )
218     do bj = jtlo,jthi
219     do bi = itlo,ithi
220     do j = jmin,jmax
221     do i = imin,imax
222     cph sss(i,j,bi,bj) = sss(i,j,bi,bj) + tmpfld2d(i,j,bi,bj)
223     salt(i,j,1,bi,bj) = salt(i,j,1,bi,bj)
224     & + tmpfld2d(i,j,bi,bj)
225     enddo
226     enddo
227     enddo
228     enddo
229     #endif
230    
231 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
232     c-- diffkr.
233     il=ilnblnk( xx_diffkr_file )
234     write(fnamediffkr(1:80),'(2a,i10.10)')
235     & xx_diffkr_file(1:il),'.',optimcycle
236 heimbach 1.11 call active_read_xyz_loc( fnamediffkr, tmpfld3d, 1,
237 heimbach 1.3 & doglobalread, ladinit, optimcycle,
238     & mythid, xx_diffkr_dummy )
239     do bj = jtlo,jthi
240     do bi = itlo,ithi
241     do k = 1,nr
242     do j = jmin,jmax
243     do i = imin,imax
244     diffkr(i,j,k,bi,bj) = diffkr(i,j,k,bi,bj) +
245     & tmpfld3d(i,j,k,bi,bj)
246     enddo
247     enddo
248     enddo
249     enddo
250     enddo
251     #endif
252    
253     #ifdef ALLOW_KAPGM_CONTROL
254     c-- kapgm.
255     il=ilnblnk( xx_kapgm_file )
256     write(fnamekapgm(1:80),'(2a,i10.10)')
257     & xx_kapgm_file(1:il),'.',optimcycle
258 heimbach 1.11 call active_read_xyz_loc( fnamekapgm, tmpfld3d, 1,
259 heimbach 1.3 & doglobalread, ladinit, optimcycle,
260     & mythid, xx_kapgm_dummy )
261     do bj = jtlo,jthi
262     do bi = itlo,ithi
263     do k = 1,nr
264     do j = jmin,jmax
265     do i = imin,imax
266     kapgm(i,j,k,bi,bj) = kapgm(i,j,k,bi,bj) +
267     & tmpfld3d(i,j,k,bi,bj)
268     enddo
269     enddo
270     enddo
271     enddo
272     enddo
273     #endif
274    
275 heimbach 1.6 #ifdef ALLOW_EFLUXY0_CONTROL
276     c-- y-component EP-flux field.
277     il=ilnblnk( xx_efluxy_file )
278     write(fnameefluxy(1:80),'(2a,i10.10)')
279     & xx_efluxy_file(1:il),'.',optimcycle
280 heimbach 1.11 call active_read_xyz_loc( fnameefluxy, tmpfld3d, 1,
281 heimbach 1.6 & doglobalread, ladinit, optimcycle,
282     & mythid, xx_efluxy_dummy )
283    
284     do bj = jtlo,jthi
285     do bi = itlo,ithi
286     do k = 1,nr
287     do j = jmin,jmax
288     do i = imin,imax
289     EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
290     & - fac*tmpfld3d(i,j,k,bi,bj)
291     & *maskS(i,j,k,bi,bj)
292     cph EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
293     cph & - rSphere*cosFacU(J,bi,bj)
294     cph & *fac*tmpfld3d(i,j,k,bi,bj)
295     enddo
296     enddo
297     enddo
298     enddo
299     enddo
300     #endif
301    
302     #ifdef ALLOW_EFLUXP0_CONTROL
303     c-- p-component EP-flux field.
304     il=ilnblnk( xx_efluxp_file )
305     write(fnameefluxp(1:80),'(2a,i10.10)')
306     & xx_efluxp_file(1:il),'.',optimcycle
307 heimbach 1.11 call active_read_xyz_loc( fnameefluxp, tmpfld3d, 1,
308 heimbach 1.6 & doglobalread, ladinit, optimcycle,
309     & mythid, xx_efluxp_dummy )
310    
311     do bj = jtlo,jthi
312     do bi = itlo,ithi
313     do k = 1,nr
314     do j = jmin,jmax
315     do i = imin,imax
316     EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
317     & + fCori(i,j,bi,bj)
318     & *fac*tmpfld3d(i,j,k,bi,bj)
319     & *hFacV(i,j,k,bi,bj)
320     cph EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
321     cph & + fCori(i,j,bi,bj)
322     cph & *rSphere*cosFacU(J,bi,bj)
323     cph & *fac*tmpfld3d(i,j,k,bi,bj)
324     enddo
325     enddo
326     enddo
327     enddo
328     enddo
329     #endif
330    
331 heimbach 1.7 #ifdef ALLOW_BOTTOMDRAG_CONTROL
332     c-- bottom drag
333     il=ilnblnk( xx_bottomdrag_file )
334     write(fnamebottomdrag(1:80),'(2a,i10.10)')
335     & xx_bottomdrag_file(1:il),'.',optimcycle
336 heimbach 1.11 call active_read_xy_loc ( fnamebottomdrag, tmpfld2d, 1,
337 heimbach 1.7 & doglobalread, ladinit, optimcycle,
338     & mythid, xx_bottomdrag_dummy )
339     do bj = jtlo,jthi
340     do bi = itlo,ithi
341     do j = jmin,jmax
342     do i = imin,imax
343     bottomdragfld(i,j,bi,bj) = bottomdragfld(i,j,bi,bj)
344     & + tmpfld2d(i,j,bi,bj)
345     enddo
346     enddo
347     enddo
348     enddo
349     #endif
350    
351 heimbach 1.15 fdef ALLOW_EDTAUX_CONTROL
352     c-- zonal eddy stress : edtaux
353     il=ilnblnk( xx_edtaux_file )
354     write(fnameedtaux(1:80),'(2a,i10.10)')
355     & xx_edtaux_file(1:il),'.',optimcycle
356     call active_read_xyz( fnameedtaux, tmpfld3d, 1,
357     & doglobalread, ladinit, optimcycle,
358     & mythid, xx_edtaux_dummy )
359     do bj = jtlo,jthi
360     do bi = itlo,ithi
361     do k = 1,nr
362     do j = jmin,jmax
363     do i = imin,imax
364     Eddytaux(i,j,k,bi,bj) = Eddytaux(i,j,k,bi,bj) +
365     & tmpfld3d(i,j,k,bi,bj)
366     enddo
367     enddo
368     enddo
369     enddo
370     enddo
371     #endif
372    
373     #ifdef ALLOW_EDTAUY_CONTROL
374     c-- meridional eddy stress : edtauy
375     il=ilnblnk( xx_edtauy_file )
376     write(fnameedtauy(1:80),'(2a,i10.10)')
377     & xx_edtauy_file(1:il),'.',optimcycle
378     call active_read_xyz( fnameedtauy, tmpfld3d, 1,
379     & doglobalread, ladinit, optimcycle,
380     & mythid, xx_edtauy_dummy )
381     do bj = jtlo,jthi
382     do bi = itlo,ithi
383     do k = 1,nr
384     do j = jmin,jmax
385     do i = imin,imax
386     Eddytauy(i,j,k,bi,bj) = Eddytauy(i,j,k,bi,bj) +
387     & tmpfld3d(i,j,k,bi,bj)
388     enddo
389     enddo
390     enddo
391     enddo
392     enddo
393     #endif
394 heimbach 1.1
395     c-- Update the tile edges.
396    
397 heimbach 1.14 #if (defined (ALLOW_THETA0_CONTROL) || defined (ALLOW_SST0_CONTROL))
398 heimbach 1.1 _EXCH_XYZ_R8( theta, mythid )
399     #endif
400 heimbach 1.14 #if (defined (ALLOW_SALT0_CONTROL) || defined (ALLOW_SSS0_CONTROL))
401 heimbach 1.1 _EXCH_XYZ_R8( salt, mythid )
402 heimbach 1.2 #endif
403     #ifdef ALLOW_TR10_CONTROL
404 heimbach 1.13 #ifdef ALLOW_PTRACERS
405     _EXCH_XYZ_R8(pTracer(1-Olx,1-Oly,1,1,1,1),myThid)
406     #endif
407 heimbach 1.1 #endif
408 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
409     _EXCH_XYZ_R8( diffkr, mythid)
410     #endif
411     #ifdef ALLOW_KAPGM_CONTROL
412     _EXCH_XYZ_R8( kapgm, mythid)
413 heimbach 1.6 #endif
414     #ifdef ALLOW_EFLUXY0_CONTROL
415     _EXCH_XYZ_R8( EfluxY, mythid )
416     #endif
417     #ifdef ALLOW_EFLUXP0_CONTROL
418     _EXCH_XYZ_R8( EfluxP, mythid )
419 heimbach 1.7 #endif
420     #ifdef ALLOW_BOTTOMDRAG_CONTROL
421     _EXCH_XY_R8( bottomdragfld, mythid )
422 heimbach 1.3 #endif
423    
424 heimbach 1.15 #if (defined (ALLOW_EDTAUX_CONTROL) && defined (ALLOW_EDTAUY_CONTROL))
425     CALL EXCH_UV_XYZ_RS(Eddytaux,Eddytauy,.TRUE.,myThid)
426     #elif (defined (ALLOW_EDTAUX_CONTROL) || defined (ALLOW_EDTAUY_CONTROL))
427     STOP 'ctrl_map_forcing: need BOTH ALLOW_EDTAU[X,Y]_CONTROL'
428     #endif
429 heimbach 1.1
430     return
431     end
432    

  ViewVC Help
Powered by ViewVC 1.1.22