/[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.45 - (hide annotations) (download)
Wed Feb 18 12:31:10 2015 UTC (10 years, 4 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65n, checkpoint65l, checkpoint65m, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65o, HEAD
Changes since 1.44: +1 -19 lines
o change to local arrays
o remove special treatment for OpenAD

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

  ViewVC Help
Powered by ViewVC 1.1.22