/[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.30 - (hide annotations) (download)
Fri May 30 02:48:28 2008 UTC (17 years, 1 month ago) by gforget
Branch: MAIN
CVS Tags: checkpoint60, checkpoint61, checkpoint61e, checkpoint61d, checkpoint61b, checkpoint61c, checkpoint61a
Changes since 1.29: +10 -15 lines
o bridging the gap between eddy stress and GM.
  -> eddyTau is replaced with eddyPsi (eddyTau = f x rho0 x eddyPsi)
      along with a change in CPP option (now ALLOW_EDDYPSI).
  -> when using GM w/ GM_AdvForm:
      The total eddy streamfunction (Psi = eddyPsi + K x Slope)
      is applied either in the tracer Eq. or in momentum Eq.
      depending on data.gmredi (intro. GM_InMomAsStress).
  -> ALLOW_EDDYPSI_CONTROL for estimation purpose.
  The key modifications are in model/src/taueddy_external_forcing.F
  pkg/gmredi/gmredi_calc_*F pkg/gmredi/gmredi_*transport.F

1 gforget 1.30 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_map_ini.F,v 1.29 2008/03/28 18:48:05 heimbach Exp $
2 jmc 1.25 C $Name: $
3 heimbach 1.1
4     #include "CTRL_CPPOPTIONS.h"
5    
6 heimbach 1.5 CBOP
7     C !ROUTINE: ctrl_map_ini
8     C !INTERFACE:
9     subroutine ctrl_map_ini( mythid )
10    
11     C !DESCRIPTION: \bv
12     c *=================================================================
13     c | SUBROUTINE ctrl_map_ini
14 jmc 1.25 c | Add the temperature, salinity, and diffusivity parts of the
15     c | control vector to the model state and update the tile halos.
16 heimbach 1.5 c | The control vector is defined in the header file "ctrl.h".
17     c *=================================================================
18     C \ev
19 heimbach 1.1
20 heimbach 1.5 C !USES:
21 heimbach 1.1 implicit none
22    
23     c == global variables ==
24 heimbach 1.6 #include "SIZE.h"
25 heimbach 1.1 #include "EEPARAMS.h"
26 heimbach 1.6 #include "PARAMS.h"
27 heimbach 1.15 #include "GRID.h"
28 heimbach 1.1 #include "DYNVARS.h"
29 heimbach 1.15 #include "FFIELDS.h"
30 heimbach 1.1 #include "ctrl.h"
31     #include "ctrl_dummy.h"
32 heimbach 1.2 #include "optim.h"
33 heimbach 1.13 #ifdef ALLOW_PTRACERS
34     # include "PTRACERS_SIZE.h"
35 jmc 1.27 c#include "PTRACERS_PARAMS.h"
36     # include "PTRACERS_FIELDS.h"
37 heimbach 1.13 #endif
38 heimbach 1.12 #ifdef ALLOW_ECCO
39     # include "ecco_cost.h"
40     #endif
41 heimbach 1.1
42 heimbach 1.5 C !INPUT/OUTPUT PARAMETERS:
43 heimbach 1.1 c == routine arguments ==
44     integer mythid
45    
46 heimbach 1.5 C !LOCAL VARIABLES:
47 heimbach 1.1 c == local variables ==
48    
49     integer bi,bj
50     integer i,j,k
51     integer itlo,ithi
52     integer jtlo,jthi
53     integer jmin,jmax
54     integer imin,imax
55     integer il
56    
57     logical equal
58     logical doglobalread
59     logical ladinit
60    
61 heimbach 1.17 character*( 80) fnamegeneric
62 heimbach 1.1
63 heimbach 1.5 _RL fac
64 heimbach 1.12 _RL tmptest
65 heimbach 1.5
66 heimbach 1.1 c == external ==
67     integer ilnblnk
68     external ilnblnk
69    
70     c == end of interface ==
71 heimbach 1.5 CEOP
72 heimbach 1.1
73     jtlo = mybylo(mythid)
74     jthi = mybyhi(mythid)
75     itlo = mybxlo(mythid)
76     ithi = mybxhi(mythid)
77 heimbach 1.8 jmin = 1
78     jmax = sny
79     imin = 1
80     imax = snx
81 heimbach 1.1
82     doglobalread = .false.
83     ladinit = .false.
84    
85     equal = .true.
86    
87     if ( equal ) then
88     fac = 1. _d 0
89     else
90     fac = 0. _d 0
91     endif
92    
93     #ifdef ALLOW_THETA0_CONTROL
94     c-- Temperature field.
95     il=ilnblnk( xx_theta_file )
96 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
97 heimbach 1.1 & xx_theta_file(1:il),'.',optimcycle
98 heimbach 1.22 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
99 heimbach 1.1 & doglobalread, ladinit, optimcycle,
100     & mythid, xx_theta_dummy )
101    
102     do bj = jtlo,jthi
103     do bi = itlo,ithi
104     do k = 1,nr
105     do j = jmin,jmax
106     do i = imin,imax
107 heimbach 1.12 #ifdef ALLOW_ECCO
108     IF (abs(tmpfld3d(i,j,k,bi,bj)).gt.
109     $ 2.0/sqrt(wtheta(k,bi,bj)))
110     $ tmpfld3d(i,j,k,bi,bj)=
111     $ sign(2.0/sqrt(wtheta(k,bi,bj)),tmpfld3d(i,j,k,bi,bj))
112     #endif
113 utke 1.23 #ifdef ALLOW_AUTODIFF_OPENAD
114 heimbach 1.18 theta(i,j,k,bi,bj) = theta(i,j,k,bi,bj) +
115 utke 1.21 & fac*xx_theta(i,j,k,bi,bj) +
116     & fac*tmpfld3d(i,j,k,bi,bj)
117 heimbach 1.18 #else
118 heimbach 1.1 theta(i,j,k,bi,bj) = theta(i,j,k,bi,bj) +
119     & fac*tmpfld3d(i,j,k,bi,bj)
120 heimbach 1.18 #endif
121 heimbach 1.26 #ifndef DISABLE_CTRL_THETA_LIMIT
122 jmc 1.25 if(theta(i,j,k,bi,bj).lt.-2.0)
123     & theta(i,j,k,bi,bj)= -2.0
124 heimbach 1.26 #endif
125 heimbach 1.1 enddo
126     enddo
127     enddo
128     enddo
129     enddo
130 heimbach 1.12
131 heimbach 1.1 #endif
132    
133     #ifdef ALLOW_SALT0_CONTROL
134     c-- Temperature field.
135     il=ilnblnk( xx_salt_file )
136 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
137 heimbach 1.1 & xx_salt_file(1:il),'.',optimcycle
138 heimbach 1.22 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
139 heimbach 1.1 & doglobalread, ladinit, optimcycle,
140     & mythid, xx_salt_dummy )
141    
142     do bj = jtlo,jthi
143     do bi = itlo,ithi
144     do k = 1,nr
145     do j = jmin,jmax
146     do i = imin,imax
147 heimbach 1.12 #ifdef ALLOW_ECCO
148     IF (abs(tmpfld3d(i,j,k,bi,bj)).gt.
149     $ 2.0/sqrt(wsalt(k,bi,bj)))
150     $ tmpfld3d(i,j,k,bi,bj)=
151     $ sign(2.0/sqrt(wsalt(k,bi,bj)),tmpfld3d(i,j,k,bi,bj))
152     #endif
153 utke 1.23 #ifdef ALLOW_AUTODIFF_OPENAD
154 heimbach 1.18 salt(i,j,k,bi,bj) = salt(i,j,k,bi,bj) +
155 utke 1.21 & fac*xx_salt(i,j,k,bi,bj) +
156     & fac*tmpfld3d(i,j,k,bi,bj)
157 heimbach 1.18 #else
158 heimbach 1.1 salt(i,j,k,bi,bj) = salt(i,j,k,bi,bj) +
159     & fac*tmpfld3d(i,j,k,bi,bj)
160 heimbach 1.18 #endif
161 heimbach 1.12
162 heimbach 1.1 enddo
163     enddo
164     enddo
165     enddo
166 heimbach 1.12 enddO
167 heimbach 1.1 #endif
168    
169 heimbach 1.2 #ifdef ALLOW_TR10_CONTROL
170 heimbach 1.13 #ifdef ALLOW_PTRACERS
171 heimbach 1.2 c-- Temperature field.
172     il=ilnblnk( xx_tr1_file )
173 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
174 heimbach 1.2 & xx_tr1_file(1:il),'.',optimcycle
175 heimbach 1.22 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
176 heimbach 1.2 & doglobalread, ladinit, optimcycle,
177     & mythid, xx_tr1_dummy )
178    
179     do bj = jtlo,jthi
180     do bi = itlo,ithi
181     do k = 1,nr
182     do j = jmin,jmax
183     do i = imin,imax
184 heimbach 1.13 ptracer(i,j,k,bi,bj,1) = ptracer(i,j,k,bi,bj,1) +
185 heimbach 1.2 & fac*tmpfld3d(i,j,k,bi,bj)
186     enddo
187     enddo
188     enddo
189     enddo
190     enddo
191     #endif
192 heimbach 1.13 #endif
193 heimbach 1.2
194 heimbach 1.14 #ifdef ALLOW_SST0_CONTROL
195     c-- sst0.
196     il=ilnblnk( xx_sst_file )
197 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
198 heimbach 1.14 & xx_sst_file(1:il),'.',optimcycle
199 heimbach 1.22 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
200 heimbach 1.14 & doglobalread, ladinit, optimcycle,
201     & mythid, xx_sst_dummy )
202     do bj = jtlo,jthi
203     do bi = itlo,ithi
204     do j = jmin,jmax
205     do i = imin,imax
206     cph sst(i,j,bi,bj) = sst(i,j,bi,bj) + tmpfld2d(i,j,bi,bj)
207 jmc 1.25 theta(i,j,1,bi,bj) = theta(i,j,1,bi,bj)
208 heimbach 1.14 & + tmpfld2d(i,j,bi,bj)
209     enddo
210     enddo
211     enddo
212     enddo
213     #endif
214    
215     #ifdef ALLOW_SSS0_CONTROL
216     c-- sss0.
217     il=ilnblnk( xx_sss_file )
218 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
219 heimbach 1.14 & xx_sss_file(1:il),'.',optimcycle
220 heimbach 1.22 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
221 heimbach 1.14 & doglobalread, ladinit, optimcycle,
222     & mythid, xx_sss_dummy )
223     do bj = jtlo,jthi
224     do bi = itlo,ithi
225     do j = jmin,jmax
226     do i = imin,imax
227     cph sss(i,j,bi,bj) = sss(i,j,bi,bj) + tmpfld2d(i,j,bi,bj)
228     salt(i,j,1,bi,bj) = salt(i,j,1,bi,bj)
229     & + tmpfld2d(i,j,bi,bj)
230     enddo
231     enddo
232     enddo
233     enddo
234     #endif
235    
236 heimbach 1.29 #ifdef ALLOW_AUTODIFF
237 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
238     c-- diffkr.
239     il=ilnblnk( xx_diffkr_file )
240 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
241 heimbach 1.3 & xx_diffkr_file(1:il),'.',optimcycle
242 heimbach 1.22 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
243 heimbach 1.3 & doglobalread, ladinit, optimcycle,
244     & mythid, xx_diffkr_dummy )
245     do bj = jtlo,jthi
246     do bi = itlo,ithi
247     do k = 1,nr
248     do j = jmin,jmax
249     do i = imin,imax
250     diffkr(i,j,k,bi,bj) = diffkr(i,j,k,bi,bj) +
251     & tmpfld3d(i,j,k,bi,bj)
252     enddo
253     enddo
254     enddo
255     enddo
256     enddo
257     #endif
258 heimbach 1.29 #endif
259 heimbach 1.3
260 heimbach 1.29 #ifdef ALLOW_AUTODIFF
261 heimbach 1.3 #ifdef ALLOW_KAPGM_CONTROL
262     c-- kapgm.
263     il=ilnblnk( xx_kapgm_file )
264 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
265 heimbach 1.3 & xx_kapgm_file(1:il),'.',optimcycle
266 heimbach 1.22 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
267 heimbach 1.3 & doglobalread, ladinit, optimcycle,
268     & mythid, xx_kapgm_dummy )
269     do bj = jtlo,jthi
270     do bi = itlo,ithi
271     do k = 1,nr
272     do j = jmin,jmax
273     do i = imin,imax
274     kapgm(i,j,k,bi,bj) = kapgm(i,j,k,bi,bj) +
275     & tmpfld3d(i,j,k,bi,bj)
276     enddo
277     enddo
278     enddo
279     enddo
280     enddo
281     #endif
282 heimbach 1.29 #endif
283 heimbach 1.3
284 heimbach 1.29 #ifdef ALLOW_AUTODIFF
285 gforget 1.28 #ifdef ALLOW_KAPREDI_CONTROL
286     c-- kapredi.
287     il=ilnblnk( xx_kapredi_file )
288     write(fnamegeneric(1:80),'(2a,i10.10)')
289     & xx_kapredi_file(1:il),'.',optimcycle
290     call active_read_xyz( fnamegeneric, tmpfld3d, 1,
291     & doglobalread, ladinit, optimcycle,
292     & mythid, xx_kapredi_dummy )
293     do bj = jtlo,jthi
294     do bi = itlo,ithi
295     do k = 1,nr
296     do j = jmin,jmax
297     do i = imin,imax
298     kapredi(i,j,k,bi,bj) = kapredi(i,j,k,bi,bj) +
299     & tmpfld3d(i,j,k,bi,bj)
300     enddo
301     enddo
302     enddo
303     enddo
304     enddo
305     #endif
306 heimbach 1.29 #endif
307 gforget 1.28
308 heimbach 1.6 #ifdef ALLOW_EFLUXY0_CONTROL
309     c-- y-component EP-flux field.
310     il=ilnblnk( xx_efluxy_file )
311 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
312 heimbach 1.6 & xx_efluxy_file(1:il),'.',optimcycle
313 heimbach 1.22 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
314 heimbach 1.6 & doglobalread, ladinit, optimcycle,
315     & mythid, xx_efluxy_dummy )
316    
317     do bj = jtlo,jthi
318     do bi = itlo,ithi
319     do k = 1,nr
320     do j = jmin,jmax
321     do i = imin,imax
322     EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
323     & - fac*tmpfld3d(i,j,k,bi,bj)
324     & *maskS(i,j,k,bi,bj)
325     cph EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
326     cph & - rSphere*cosFacU(J,bi,bj)
327     cph & *fac*tmpfld3d(i,j,k,bi,bj)
328     enddo
329     enddo
330     enddo
331     enddo
332     enddo
333     #endif
334    
335     #ifdef ALLOW_EFLUXP0_CONTROL
336     c-- p-component EP-flux field.
337     il=ilnblnk( xx_efluxp_file )
338 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
339 heimbach 1.6 & xx_efluxp_file(1:il),'.',optimcycle
340 heimbach 1.22 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
341 heimbach 1.6 & doglobalread, ladinit, optimcycle,
342     & mythid, xx_efluxp_dummy )
343    
344     do bj = jtlo,jthi
345     do bi = itlo,ithi
346     do k = 1,nr
347     do j = jmin,jmax
348     do i = imin,imax
349     EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
350     & + fCori(i,j,bi,bj)
351     & *fac*tmpfld3d(i,j,k,bi,bj)
352     & *hFacV(i,j,k,bi,bj)
353     cph EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
354     cph & + fCori(i,j,bi,bj)
355     cph & *rSphere*cosFacU(J,bi,bj)
356     cph & *fac*tmpfld3d(i,j,k,bi,bj)
357     enddo
358     enddo
359     enddo
360     enddo
361     enddo
362     #endif
363    
364 heimbach 1.7 #ifdef ALLOW_BOTTOMDRAG_CONTROL
365     c-- bottom drag
366     il=ilnblnk( xx_bottomdrag_file )
367 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
368 heimbach 1.7 & xx_bottomdrag_file(1:il),'.',optimcycle
369 heimbach 1.22 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
370 heimbach 1.7 & doglobalread, ladinit, optimcycle,
371     & mythid, xx_bottomdrag_dummy )
372     do bj = jtlo,jthi
373     do bi = itlo,ithi
374     do j = jmin,jmax
375     do i = imin,imax
376 jmc 1.25 bottomdragfld(i,j,bi,bj) = bottomdragfld(i,j,bi,bj)
377 heimbach 1.7 & + tmpfld2d(i,j,bi,bj)
378     enddo
379     enddo
380     enddo
381     enddo
382     #endif
383    
384 gforget 1.30 #ifdef ALLOW_EDDYPSI_CONTROL
385     c-- zonal eddy streamfunction : eddyPsiX
386 heimbach 1.15 il=ilnblnk( xx_edtaux_file )
387 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
388 heimbach 1.15 & xx_edtaux_file(1:il),'.',optimcycle
389 heimbach 1.17 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
390 heimbach 1.15 & doglobalread, ladinit, optimcycle,
391     & mythid, xx_edtaux_dummy )
392     do bj = jtlo,jthi
393     do bi = itlo,ithi
394     do k = 1,nr
395     do j = jmin,jmax
396     do i = imin,imax
397 gforget 1.30 eddyPsiX(i,j,k,bi,bj) = eddyPsiX(i,j,k,bi,bj) +
398     & tmpfld3d(i,j,k,bi,bj)
399 heimbach 1.15 enddo
400     enddo
401     enddo
402     enddo
403     enddo
404 gforget 1.30 c-- meridional eddy streamfunction : eddyPsiY
405 heimbach 1.15 il=ilnblnk( xx_edtauy_file )
406 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
407 heimbach 1.15 & xx_edtauy_file(1:il),'.',optimcycle
408 heimbach 1.17 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
409 heimbach 1.15 & doglobalread, ladinit, optimcycle,
410     & mythid, xx_edtauy_dummy )
411     do bj = jtlo,jthi
412     do bi = itlo,ithi
413     do k = 1,nr
414     do j = jmin,jmax
415     do i = imin,imax
416 gforget 1.30 eddyPsiY(i,j,k,bi,bj) = eddyPsiY(i,j,k,bi,bj) +
417     & tmpfld3d(i,j,k,bi,bj)
418 heimbach 1.15 enddo
419     enddo
420     enddo
421     enddo
422     enddo
423     #endif
424 heimbach 1.1
425 heimbach 1.17 #ifdef ALLOW_UVEL0_CONTROL
426     c-- initial zonal velocity
427     il=ilnblnk( xx_uvel_file )
428     write(fnamegeneric(1:80),'(2a,i10.10)')
429     & xx_uvel_file(1:il),'.',optimcycle
430     call active_read_xyz( fnamegeneric, tmpfld3d, 1,
431     & doglobalread, ladinit, optimcycle,
432     & mythid, xx_uvel_dummy )
433     do bj = jtlo,jthi
434     do bi = itlo,ithi
435     do k = 1,nr
436     do j = jmin,jmax
437     do i = imin,imax
438 utke 1.23 #ifdef ALLOW_AUTODIFF_OPENAD
439 heimbach 1.19 uVel(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) +
440     & fac*xx_uvel(i,j,k,bi,bj)
441     #else
442 heimbach 1.17 uVel(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) +
443 heimbach 1.19 & fac*tmpfld3d(i,j,k,bi,bj)
444     #endif
445 heimbach 1.17 enddo
446     enddo
447     enddo
448     enddo
449     enddo
450     #endif
451    
452     #ifdef ALLOW_VVEL0_CONTROL
453     c-- initial merid. velocity
454     il=ilnblnk( xx_vvel_file )
455     write(fnamegeneric(1:80),'(2a,i10.10)')
456     & xx_vvel_file(1:il),'.',optimcycle
457     call active_read_xyz( fnamegeneric, tmpfld3d, 1,
458     & doglobalread, ladinit, optimcycle,
459     & mythid, xx_vvel_dummy )
460     do bj = jtlo,jthi
461     do bi = itlo,ithi
462     do k = 1,nr
463     do j = jmin,jmax
464     do i = imin,imax
465 utke 1.23 #ifdef ALLOW_AUTODIFF_OPENAD
466 heimbach 1.17 vVel(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) +
467 heimbach 1.19 & fac*xx_vvel(i,j,k,bi,bj)
468     #else
469     vVel(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) +
470     & fac*tmpfld3d(i,j,k,bi,bj)
471     #endif
472 heimbach 1.17 enddo
473     enddo
474     enddo
475     enddo
476     enddo
477     #endif
478    
479     #ifdef ALLOW_ETAN0_CONTROL
480     c-- initial Eta.
481     il=ilnblnk( xx_etan_file )
482     write(fnamegeneric(1:80),'(2a,i10.10)')
483     & xx_etan_file(1:il),'.',optimcycle
484 heimbach 1.22 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
485 heimbach 1.17 & doglobalread, ladinit, optimcycle,
486     & mythid, xx_etan_dummy )
487     do bj = jtlo,jthi
488     do bi = itlo,ithi
489     do j = jmin,jmax
490     do i = imin,imax
491 utke 1.23 #ifdef ALLOW_AUTODIFF_OPENAD
492 heimbach 1.19 etaN(i,j,bi,bj) = etaN(i,j,bi,bj) +
493     & fac*xx_etan(i,j,bi,bj)
494     #else
495 jmc 1.25 etaN(i,j,bi,bj) = etaN(i,j,bi,bj) +
496 heimbach 1.19 & fac*tmpfld2d(i,j,bi,bj)
497     #endif
498 heimbach 1.17 enddo
499     enddo
500     enddo
501     enddo
502     #endif
503    
504     #ifdef ALLOW_RELAXSST_CONTROL
505     c-- SST relaxation coefficient.
506     il=ilnblnk( xx_relaxsst_file )
507     write(fnamegeneric(1:80),'(2a,i10.10)')
508     & xx_relaxsst_file(1:il),'.',optimcycle
509 heimbach 1.22 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
510 heimbach 1.17 & doglobalread, ladinit, optimcycle,
511     & mythid, xx_relaxsst_dummy )
512     do bj = jtlo,jthi
513     do bi = itlo,ithi
514     do j = jmin,jmax
515     do i = imin,imax
516 jmc 1.25 lambdaThetaClimRelax(i,j,bi,bj) =
517     & lambdaThetaClimRelax(i,j,bi,bj)
518 heimbach 1.17 & + tmpfld2d(i,j,bi,bj)
519     enddo
520     enddo
521     enddo
522     enddo
523     #endif
524    
525     #ifdef ALLOW_RELAXSSS_CONTROL
526     c-- SSS relaxation coefficient.
527     il=ilnblnk( xx_relaxsss_file )
528     write(fnamegeneric(1:80),'(2a,i10.10)')
529     & xx_relaxsss_file(1:il),'.',optimcycle
530 heimbach 1.22 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
531 heimbach 1.17 & doglobalread, ladinit, optimcycle,
532     & mythid, xx_relaxsss_dummy )
533     do bj = jtlo,jthi
534     do bi = itlo,ithi
535     do j = jmin,jmax
536     do i = imin,imax
537 jmc 1.25 lambdaSaltClimRelax(i,j,bi,bj) =
538     & lambdaSaltClimRelax(i,j,bi,bj)
539 heimbach 1.17 & + tmpfld2d(i,j,bi,bj)
540     enddo
541     enddo
542     enddo
543     enddo
544     #endif
545    
546 heimbach 1.24 #ifdef ALLOW_SEAICE
547     call seaice_ctrl_map_ini( mythid )
548     #endif
549    
550 heimbach 1.1 c-- Update the tile edges.
551    
552 heimbach 1.14 #if (defined (ALLOW_THETA0_CONTROL) || defined (ALLOW_SST0_CONTROL))
553 heimbach 1.1 _EXCH_XYZ_R8( theta, mythid )
554     #endif
555 heimbach 1.14 #if (defined (ALLOW_SALT0_CONTROL) || defined (ALLOW_SSS0_CONTROL))
556 heimbach 1.1 _EXCH_XYZ_R8( salt, mythid )
557 heimbach 1.2 #endif
558     #ifdef ALLOW_TR10_CONTROL
559 heimbach 1.13 #ifdef ALLOW_PTRACERS
560     _EXCH_XYZ_R8(pTracer(1-Olx,1-Oly,1,1,1,1),myThid)
561     #endif
562 heimbach 1.1 #endif
563 heimbach 1.29
564     #ifdef ALLOW_AUTODIFF
565     # ifdef ALLOW_DIFFKR_CONTROL
566 heimbach 1.3 _EXCH_XYZ_R8( diffkr, mythid)
567 heimbach 1.29 # endif
568     # ifdef ALLOW_KAPGM_CONTROL
569 heimbach 1.3 _EXCH_XYZ_R8( kapgm, mythid)
570 heimbach 1.29 # endif
571     # ifdef ALLOW_KAPREDI_CONTROL
572 gforget 1.28 _EXCH_XYZ_R8( kapredi, mythid)
573 heimbach 1.29 # endif
574 gforget 1.28 #endif
575 heimbach 1.29
576 heimbach 1.6 #ifdef ALLOW_EFLUXY0_CONTROL
577     _EXCH_XYZ_R8( EfluxY, mythid )
578     #endif
579     #ifdef ALLOW_EFLUXP0_CONTROL
580     _EXCH_XYZ_R8( EfluxP, mythid )
581 heimbach 1.7 #endif
582     #ifdef ALLOW_BOTTOMDRAG_CONTROL
583     _EXCH_XY_R8( bottomdragfld, mythid )
584 heimbach 1.3 #endif
585    
586 gforget 1.30 #ifdef ALLOW_EDDYPSI_CONTROL
587     CALL EXCH_UV_XYZ_RS(eddyPsiX,eddyPsiY,.TRUE.,myThid)
588 heimbach 1.15 #endif
589 heimbach 1.1
590 heimbach 1.17 #ifdef ALLOW_UVEL0_CONTROL
591     _EXCH_XYZ_R8( uVel, mythid)
592     #endif
593    
594     #ifdef ALLOW_VVEL0_CONTROL
595     _EXCH_XYZ_R8( vVel, mythid)
596     #endif
597    
598     #ifdef ALLOW_ETAN0_CONTROL
599     _EXCH_XY_R8( etaN, mythid )
600     #endif
601    
602     #ifdef ALLOW_RELAXSST_CONTROL
603     _EXCH_XY_R4( lambdaThetaClimRelax, mythid )
604     #endif
605    
606     #ifdef ALLOW_RELAXSSS_CONTROL
607     _EXCH_XY_R4( lambdaThetaClimRelax, mythid )
608     #endif
609    
610 heimbach 1.1 return
611     end
612    

  ViewVC Help
Powered by ViewVC 1.1.22