/[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.32 - (hide annotations) (download)
Tue Apr 28 18:09:28 2009 UTC (15 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62t, checkpoint62c, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62, checkpoint62b, checkpoint61n, checkpoint61q, checkpoint61o, checkpoint61m, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.31: +15 -15 lines
change macros (EXCH & GLOBAL_SUM/MAX) sufix _R4/_R8 to _RS/_RL
 when applied to _RS/_RL variable

1 jmc 1.32 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_map_ini.F,v 1.31 2008/11/18 16:47:15 utke 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 utke 1.31 #ifdef ALLOW_AUTODIFF_OPENAD
251 heimbach 1.3 diffkr(i,j,k,bi,bj) = diffkr(i,j,k,bi,bj) +
252 utke 1.31 & xx_diffkr(i,j,k,bi,bj) +
253 heimbach 1.3 & tmpfld3d(i,j,k,bi,bj)
254 utke 1.31 #else
255     diffkr(i,j,k,bi,bj) = diffkr(i,j,k,bi,bj) +
256     & tmpfld3d(i,j,k,bi,bj)
257     #endif
258 heimbach 1.3 enddo
259     enddo
260     enddo
261     enddo
262     enddo
263     #endif
264 heimbach 1.29 #endif
265 heimbach 1.3
266 heimbach 1.29 #ifdef ALLOW_AUTODIFF
267 heimbach 1.3 #ifdef ALLOW_KAPGM_CONTROL
268     c-- kapgm.
269     il=ilnblnk( xx_kapgm_file )
270 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
271 heimbach 1.3 & xx_kapgm_file(1:il),'.',optimcycle
272 heimbach 1.22 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
273 heimbach 1.3 & doglobalread, ladinit, optimcycle,
274     & mythid, xx_kapgm_dummy )
275     do bj = jtlo,jthi
276     do bi = itlo,ithi
277     do k = 1,nr
278     do j = jmin,jmax
279     do i = imin,imax
280 utke 1.31 #ifdef ALLOW_AUTODIFF_OPENAD
281 heimbach 1.3 kapgm(i,j,k,bi,bj) = kapgm(i,j,k,bi,bj) +
282 utke 1.31 & xx_kapgm(i,j,k,bi,bj) +
283 heimbach 1.3 & tmpfld3d(i,j,k,bi,bj)
284 utke 1.31 #else
285     kapgm(i,j,k,bi,bj) = kapgm(i,j,k,bi,bj) +
286     & tmpfld3d(i,j,k,bi,bj)
287     #endif
288 heimbach 1.3 enddo
289     enddo
290     enddo
291     enddo
292     enddo
293     #endif
294 heimbach 1.29 #endif
295 heimbach 1.3
296 heimbach 1.29 #ifdef ALLOW_AUTODIFF
297 gforget 1.28 #ifdef ALLOW_KAPREDI_CONTROL
298     c-- kapredi.
299     il=ilnblnk( xx_kapredi_file )
300     write(fnamegeneric(1:80),'(2a,i10.10)')
301     & xx_kapredi_file(1:il),'.',optimcycle
302     call active_read_xyz( fnamegeneric, tmpfld3d, 1,
303     & doglobalread, ladinit, optimcycle,
304     & mythid, xx_kapredi_dummy )
305     do bj = jtlo,jthi
306     do bi = itlo,ithi
307     do k = 1,nr
308     do j = jmin,jmax
309     do i = imin,imax
310     kapredi(i,j,k,bi,bj) = kapredi(i,j,k,bi,bj) +
311     & tmpfld3d(i,j,k,bi,bj)
312     enddo
313     enddo
314     enddo
315     enddo
316     enddo
317     #endif
318 heimbach 1.29 #endif
319 gforget 1.28
320 heimbach 1.6 #ifdef ALLOW_EFLUXY0_CONTROL
321     c-- y-component EP-flux field.
322     il=ilnblnk( xx_efluxy_file )
323 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
324 heimbach 1.6 & xx_efluxy_file(1:il),'.',optimcycle
325 heimbach 1.22 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
326 heimbach 1.6 & doglobalread, ladinit, optimcycle,
327     & mythid, xx_efluxy_dummy )
328    
329     do bj = jtlo,jthi
330     do bi = itlo,ithi
331     do k = 1,nr
332     do j = jmin,jmax
333     do i = imin,imax
334     EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
335     & - fac*tmpfld3d(i,j,k,bi,bj)
336     & *maskS(i,j,k,bi,bj)
337     cph EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
338     cph & - rSphere*cosFacU(J,bi,bj)
339     cph & *fac*tmpfld3d(i,j,k,bi,bj)
340     enddo
341     enddo
342     enddo
343     enddo
344     enddo
345     #endif
346    
347     #ifdef ALLOW_EFLUXP0_CONTROL
348     c-- p-component EP-flux field.
349     il=ilnblnk( xx_efluxp_file )
350 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
351 heimbach 1.6 & xx_efluxp_file(1:il),'.',optimcycle
352 heimbach 1.22 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
353 heimbach 1.6 & doglobalread, ladinit, optimcycle,
354     & mythid, xx_efluxp_dummy )
355    
356     do bj = jtlo,jthi
357     do bi = itlo,ithi
358     do k = 1,nr
359     do j = jmin,jmax
360     do i = imin,imax
361     EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
362     & + fCori(i,j,bi,bj)
363     & *fac*tmpfld3d(i,j,k,bi,bj)
364     & *hFacV(i,j,k,bi,bj)
365     cph EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
366     cph & + fCori(i,j,bi,bj)
367     cph & *rSphere*cosFacU(J,bi,bj)
368     cph & *fac*tmpfld3d(i,j,k,bi,bj)
369     enddo
370     enddo
371     enddo
372     enddo
373     enddo
374     #endif
375    
376 heimbach 1.7 #ifdef ALLOW_BOTTOMDRAG_CONTROL
377     c-- bottom drag
378     il=ilnblnk( xx_bottomdrag_file )
379 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
380 heimbach 1.7 & xx_bottomdrag_file(1:il),'.',optimcycle
381 heimbach 1.22 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
382 heimbach 1.7 & doglobalread, ladinit, optimcycle,
383     & mythid, xx_bottomdrag_dummy )
384     do bj = jtlo,jthi
385     do bi = itlo,ithi
386     do j = jmin,jmax
387     do i = imin,imax
388 jmc 1.25 bottomdragfld(i,j,bi,bj) = bottomdragfld(i,j,bi,bj)
389 heimbach 1.7 & + tmpfld2d(i,j,bi,bj)
390     enddo
391     enddo
392     enddo
393     enddo
394     #endif
395    
396 gforget 1.30 #ifdef ALLOW_EDDYPSI_CONTROL
397     c-- zonal eddy streamfunction : eddyPsiX
398 heimbach 1.15 il=ilnblnk( xx_edtaux_file )
399 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
400 heimbach 1.15 & xx_edtaux_file(1:il),'.',optimcycle
401 heimbach 1.17 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
402 heimbach 1.15 & doglobalread, ladinit, optimcycle,
403     & mythid, xx_edtaux_dummy )
404     do bj = jtlo,jthi
405     do bi = itlo,ithi
406     do k = 1,nr
407     do j = jmin,jmax
408     do i = imin,imax
409 gforget 1.30 eddyPsiX(i,j,k,bi,bj) = eddyPsiX(i,j,k,bi,bj) +
410     & tmpfld3d(i,j,k,bi,bj)
411 heimbach 1.15 enddo
412     enddo
413     enddo
414     enddo
415     enddo
416 gforget 1.30 c-- meridional eddy streamfunction : eddyPsiY
417 heimbach 1.15 il=ilnblnk( xx_edtauy_file )
418 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
419 heimbach 1.15 & xx_edtauy_file(1:il),'.',optimcycle
420 heimbach 1.17 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
421 heimbach 1.15 & doglobalread, ladinit, optimcycle,
422     & mythid, xx_edtauy_dummy )
423     do bj = jtlo,jthi
424     do bi = itlo,ithi
425     do k = 1,nr
426     do j = jmin,jmax
427     do i = imin,imax
428 gforget 1.30 eddyPsiY(i,j,k,bi,bj) = eddyPsiY(i,j,k,bi,bj) +
429     & tmpfld3d(i,j,k,bi,bj)
430 heimbach 1.15 enddo
431     enddo
432     enddo
433     enddo
434     enddo
435     #endif
436 heimbach 1.1
437 heimbach 1.17 #ifdef ALLOW_UVEL0_CONTROL
438     c-- initial zonal velocity
439     il=ilnblnk( xx_uvel_file )
440     write(fnamegeneric(1:80),'(2a,i10.10)')
441     & xx_uvel_file(1:il),'.',optimcycle
442     call active_read_xyz( fnamegeneric, tmpfld3d, 1,
443     & doglobalread, ladinit, optimcycle,
444     & mythid, xx_uvel_dummy )
445     do bj = jtlo,jthi
446     do bi = itlo,ithi
447     do k = 1,nr
448     do j = jmin,jmax
449     do i = imin,imax
450 utke 1.23 #ifdef ALLOW_AUTODIFF_OPENAD
451 heimbach 1.19 uVel(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) +
452     & fac*xx_uvel(i,j,k,bi,bj)
453     #else
454 heimbach 1.17 uVel(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) +
455 heimbach 1.19 & fac*tmpfld3d(i,j,k,bi,bj)
456     #endif
457 heimbach 1.17 enddo
458     enddo
459     enddo
460     enddo
461     enddo
462     #endif
463    
464     #ifdef ALLOW_VVEL0_CONTROL
465     c-- initial merid. velocity
466     il=ilnblnk( xx_vvel_file )
467     write(fnamegeneric(1:80),'(2a,i10.10)')
468     & xx_vvel_file(1:il),'.',optimcycle
469     call active_read_xyz( fnamegeneric, tmpfld3d, 1,
470     & doglobalread, ladinit, optimcycle,
471     & mythid, xx_vvel_dummy )
472     do bj = jtlo,jthi
473     do bi = itlo,ithi
474     do k = 1,nr
475     do j = jmin,jmax
476     do i = imin,imax
477 utke 1.23 #ifdef ALLOW_AUTODIFF_OPENAD
478 heimbach 1.17 vVel(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) +
479 heimbach 1.19 & fac*xx_vvel(i,j,k,bi,bj)
480     #else
481     vVel(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) +
482     & fac*tmpfld3d(i,j,k,bi,bj)
483     #endif
484 heimbach 1.17 enddo
485     enddo
486     enddo
487     enddo
488     enddo
489     #endif
490    
491     #ifdef ALLOW_ETAN0_CONTROL
492     c-- initial Eta.
493     il=ilnblnk( xx_etan_file )
494     write(fnamegeneric(1:80),'(2a,i10.10)')
495     & xx_etan_file(1:il),'.',optimcycle
496 heimbach 1.22 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
497 heimbach 1.17 & doglobalread, ladinit, optimcycle,
498     & mythid, xx_etan_dummy )
499     do bj = jtlo,jthi
500     do bi = itlo,ithi
501     do j = jmin,jmax
502     do i = imin,imax
503 utke 1.23 #ifdef ALLOW_AUTODIFF_OPENAD
504 heimbach 1.19 etaN(i,j,bi,bj) = etaN(i,j,bi,bj) +
505     & fac*xx_etan(i,j,bi,bj)
506     #else
507 jmc 1.25 etaN(i,j,bi,bj) = etaN(i,j,bi,bj) +
508 heimbach 1.19 & fac*tmpfld2d(i,j,bi,bj)
509     #endif
510 heimbach 1.17 enddo
511     enddo
512     enddo
513     enddo
514     #endif
515    
516     #ifdef ALLOW_RELAXSST_CONTROL
517     c-- SST relaxation coefficient.
518     il=ilnblnk( xx_relaxsst_file )
519     write(fnamegeneric(1:80),'(2a,i10.10)')
520     & xx_relaxsst_file(1:il),'.',optimcycle
521 heimbach 1.22 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
522 heimbach 1.17 & doglobalread, ladinit, optimcycle,
523     & mythid, xx_relaxsst_dummy )
524     do bj = jtlo,jthi
525     do bi = itlo,ithi
526     do j = jmin,jmax
527     do i = imin,imax
528 jmc 1.25 lambdaThetaClimRelax(i,j,bi,bj) =
529     & lambdaThetaClimRelax(i,j,bi,bj)
530 heimbach 1.17 & + tmpfld2d(i,j,bi,bj)
531     enddo
532     enddo
533     enddo
534     enddo
535     #endif
536    
537     #ifdef ALLOW_RELAXSSS_CONTROL
538     c-- SSS relaxation coefficient.
539     il=ilnblnk( xx_relaxsss_file )
540     write(fnamegeneric(1:80),'(2a,i10.10)')
541     & xx_relaxsss_file(1:il),'.',optimcycle
542 heimbach 1.22 call active_read_xy ( fnamegeneric, tmpfld2d, 1,
543 heimbach 1.17 & doglobalread, ladinit, optimcycle,
544     & mythid, xx_relaxsss_dummy )
545     do bj = jtlo,jthi
546     do bi = itlo,ithi
547     do j = jmin,jmax
548     do i = imin,imax
549 jmc 1.25 lambdaSaltClimRelax(i,j,bi,bj) =
550     & lambdaSaltClimRelax(i,j,bi,bj)
551 heimbach 1.17 & + tmpfld2d(i,j,bi,bj)
552     enddo
553     enddo
554     enddo
555     enddo
556     #endif
557    
558 heimbach 1.24 #ifdef ALLOW_SEAICE
559     call seaice_ctrl_map_ini( mythid )
560     #endif
561    
562 heimbach 1.1 c-- Update the tile edges.
563    
564 heimbach 1.14 #if (defined (ALLOW_THETA0_CONTROL) || defined (ALLOW_SST0_CONTROL))
565 jmc 1.32 _EXCH_XYZ_RL( theta, mythid )
566 heimbach 1.1 #endif
567 heimbach 1.14 #if (defined (ALLOW_SALT0_CONTROL) || defined (ALLOW_SSS0_CONTROL))
568 jmc 1.32 _EXCH_XYZ_RL( salt, mythid )
569 heimbach 1.2 #endif
570     #ifdef ALLOW_TR10_CONTROL
571 heimbach 1.13 #ifdef ALLOW_PTRACERS
572 jmc 1.32 _EXCH_XYZ_RL(pTracer(1-Olx,1-Oly,1,1,1,1),myThid)
573 heimbach 1.13 #endif
574 heimbach 1.1 #endif
575 heimbach 1.29
576     #ifdef ALLOW_AUTODIFF
577     # ifdef ALLOW_DIFFKR_CONTROL
578 jmc 1.32 _EXCH_XYZ_RL( diffkr, mythid)
579 heimbach 1.29 # endif
580     # ifdef ALLOW_KAPGM_CONTROL
581 jmc 1.32 _EXCH_XYZ_RL( kapgm, mythid)
582 heimbach 1.29 # endif
583     # ifdef ALLOW_KAPREDI_CONTROL
584 jmc 1.32 _EXCH_XYZ_RL( kapredi, mythid)
585 heimbach 1.29 # endif
586 gforget 1.28 #endif
587 heimbach 1.29
588 heimbach 1.6 #ifdef ALLOW_EFLUXY0_CONTROL
589 jmc 1.32 _EXCH_XYZ_RL( EfluxY, mythid )
590 heimbach 1.6 #endif
591     #ifdef ALLOW_EFLUXP0_CONTROL
592 jmc 1.32 _EXCH_XYZ_RL( EfluxP, mythid )
593 heimbach 1.7 #endif
594     #ifdef ALLOW_BOTTOMDRAG_CONTROL
595 jmc 1.32 _EXCH_XY_RL( bottomdragfld, mythid )
596 heimbach 1.3 #endif
597    
598 gforget 1.30 #ifdef ALLOW_EDDYPSI_CONTROL
599     CALL EXCH_UV_XYZ_RS(eddyPsiX,eddyPsiY,.TRUE.,myThid)
600 heimbach 1.15 #endif
601 heimbach 1.1
602 heimbach 1.17 #ifdef ALLOW_UVEL0_CONTROL
603 jmc 1.32 _EXCH_XYZ_RL( uVel, mythid)
604 heimbach 1.17 #endif
605    
606     #ifdef ALLOW_VVEL0_CONTROL
607 jmc 1.32 _EXCH_XYZ_RL( vVel, mythid)
608 heimbach 1.17 #endif
609    
610     #ifdef ALLOW_ETAN0_CONTROL
611 jmc 1.32 _EXCH_XY_RL( etaN, mythid )
612 heimbach 1.17 #endif
613    
614     #ifdef ALLOW_RELAXSST_CONTROL
615 jmc 1.32 _EXCH_XY_RS( lambdaThetaClimRelax, mythid )
616 heimbach 1.17 #endif
617    
618     #ifdef ALLOW_RELAXSSS_CONTROL
619 jmc 1.32 _EXCH_XY_RS( lambdaThetaClimRelax, mythid )
620 heimbach 1.17 #endif
621    
622 heimbach 1.1 return
623     end
624    

  ViewVC Help
Powered by ViewVC 1.1.22