/[MITgcm]/MITgcm_contrib/dgoldberg/streamice/ctrl_map_ini.F
ViewVC logotype

Annotation of /MITgcm_contrib/dgoldberg/streamice/ctrl_map_ini.F

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


Revision 1.2 - (hide annotations) (download)
Fri Jul 20 13:40:37 2012 UTC (13 years ago) by heimbach
Branch: MAIN
CVS Tags: HEAD
Changes since 1.1: +1 -1 lines
FILE REMOVED
Move some files

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

  ViewVC Help
Powered by ViewVC 1.1.22