/[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.19 - (hide annotations) (download)
Thu Jan 5 17:48:01 2006 UTC (18 years, 4 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58
Changes since 1.18: +20 -4 lines
More changes.

1 heimbach 1.19 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_map_ini.F,v 1.18 2006/01/03 17:10:35 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 heimbach 1.17 character*( 80) fnamegeneric
60 heimbach 1.1
61 heimbach 1.5 _RL fac
62 heimbach 1.12 _RL tmptest
63 heimbach 1.5
64 heimbach 1.1 c == external ==
65     integer ilnblnk
66     external ilnblnk
67    
68     c == end of interface ==
69 heimbach 1.5 CEOP
70 heimbach 1.1
71     jtlo = mybylo(mythid)
72     jthi = mybyhi(mythid)
73     itlo = mybxlo(mythid)
74     ithi = mybxhi(mythid)
75 heimbach 1.8 jmin = 1
76     jmax = sny
77     imin = 1
78     imax = snx
79 heimbach 1.1
80     doglobalread = .false.
81     ladinit = .false.
82    
83     equal = .true.
84    
85     if ( equal ) then
86     fac = 1. _d 0
87     else
88     fac = 0. _d 0
89     endif
90    
91     #ifdef ALLOW_THETA0_CONTROL
92     c-- Temperature field.
93     il=ilnblnk( xx_theta_file )
94 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
95 heimbach 1.1 & xx_theta_file(1:il),'.',optimcycle
96 heimbach 1.17 call active_read_xyz_loc( fnamegeneric, tmpfld3d, 1,
97 heimbach 1.1 & doglobalread, ladinit, optimcycle,
98     & mythid, xx_theta_dummy )
99    
100     do bj = jtlo,jthi
101     do bi = itlo,ithi
102     do k = 1,nr
103     do j = jmin,jmax
104     do i = imin,imax
105 heimbach 1.12 #ifdef ALLOW_ECCO
106     IF (abs(tmpfld3d(i,j,k,bi,bj)).gt.
107     $ 2.0/sqrt(wtheta(k,bi,bj)))
108     $ tmpfld3d(i,j,k,bi,bj)=
109     $ sign(2.0/sqrt(wtheta(k,bi,bj)),tmpfld3d(i,j,k,bi,bj))
110     #endif
111 heimbach 1.18 #ifdef ALLOW_OPENAD
112     theta(i,j,k,bi,bj) = theta(i,j,k,bi,bj) +
113     & fac*xx_theta(i,j,k,bi,bj)
114     #else
115 heimbach 1.1 theta(i,j,k,bi,bj) = theta(i,j,k,bi,bj) +
116     & fac*tmpfld3d(i,j,k,bi,bj)
117 heimbach 1.18 #endif
118 heimbach 1.9 if(theta(i,j,k,bi,bj).lt.-2.0)
119 heimbach 1.12 & theta(i,j,k,bi,bj)= -2.0
120 heimbach 1.1 enddo
121     enddo
122     enddo
123     enddo
124     enddo
125 heimbach 1.12
126 heimbach 1.1 #endif
127    
128     #ifdef ALLOW_SALT0_CONTROL
129     c-- Temperature field.
130     il=ilnblnk( xx_salt_file )
131 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
132 heimbach 1.1 & xx_salt_file(1:il),'.',optimcycle
133 heimbach 1.17 call active_read_xyz_loc( fnamegeneric, tmpfld3d, 1,
134 heimbach 1.1 & doglobalread, ladinit, optimcycle,
135     & mythid, xx_salt_dummy )
136    
137     do bj = jtlo,jthi
138     do bi = itlo,ithi
139     do k = 1,nr
140     do j = jmin,jmax
141     do i = imin,imax
142 heimbach 1.12 #ifdef ALLOW_ECCO
143     IF (abs(tmpfld3d(i,j,k,bi,bj)).gt.
144     $ 2.0/sqrt(wsalt(k,bi,bj)))
145     $ tmpfld3d(i,j,k,bi,bj)=
146     $ sign(2.0/sqrt(wsalt(k,bi,bj)),tmpfld3d(i,j,k,bi,bj))
147     #endif
148 heimbach 1.18 #ifdef ALLOW_OPENAD
149     salt(i,j,k,bi,bj) = salt(i,j,k,bi,bj) +
150     & fac*xx_salt(i,j,k,bi,bj)
151     #else
152 heimbach 1.1 salt(i,j,k,bi,bj) = salt(i,j,k,bi,bj) +
153     & fac*tmpfld3d(i,j,k,bi,bj)
154 heimbach 1.18 #endif
155 heimbach 1.12
156 heimbach 1.1 enddo
157     enddo
158     enddo
159     enddo
160 heimbach 1.12 enddO
161 heimbach 1.1 #endif
162    
163 heimbach 1.2 #ifdef ALLOW_TR10_CONTROL
164 heimbach 1.13 #ifdef ALLOW_PTRACERS
165 heimbach 1.2 c-- Temperature field.
166     il=ilnblnk( xx_tr1_file )
167 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
168 heimbach 1.2 & xx_tr1_file(1:il),'.',optimcycle
169 heimbach 1.17 call active_read_xyz_loc( fnamegeneric, tmpfld3d, 1,
170 heimbach 1.2 & doglobalread, ladinit, optimcycle,
171     & mythid, xx_tr1_dummy )
172    
173     do bj = jtlo,jthi
174     do bi = itlo,ithi
175     do k = 1,nr
176     do j = jmin,jmax
177     do i = imin,imax
178 heimbach 1.13 ptracer(i,j,k,bi,bj,1) = ptracer(i,j,k,bi,bj,1) +
179 heimbach 1.2 & fac*tmpfld3d(i,j,k,bi,bj)
180     enddo
181     enddo
182     enddo
183     enddo
184     enddo
185     #endif
186 heimbach 1.13 #endif
187 heimbach 1.2
188 heimbach 1.14 #ifdef ALLOW_SST0_CONTROL
189     c-- sst0.
190     il=ilnblnk( xx_sst_file )
191 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
192 heimbach 1.14 & xx_sst_file(1:il),'.',optimcycle
193 heimbach 1.17 call active_read_xy_loc ( fnamegeneric, tmpfld2d, 1,
194 heimbach 1.14 & doglobalread, ladinit, optimcycle,
195     & mythid, xx_sst_dummy )
196     do bj = jtlo,jthi
197     do bi = itlo,ithi
198     do j = jmin,jmax
199     do i = imin,imax
200     cph sst(i,j,bi,bj) = sst(i,j,bi,bj) + tmpfld2d(i,j,bi,bj)
201     theta(i,j,1,bi,bj) = theta(i,j,1,bi,bj)
202     & + tmpfld2d(i,j,bi,bj)
203     enddo
204     enddo
205     enddo
206     enddo
207     #endif
208    
209     #ifdef ALLOW_SSS0_CONTROL
210     c-- sss0.
211     il=ilnblnk( xx_sss_file )
212 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
213 heimbach 1.14 & xx_sss_file(1:il),'.',optimcycle
214 heimbach 1.17 call active_read_xy_loc ( fnamegeneric, tmpfld2d, 1,
215 heimbach 1.14 & doglobalread, ladinit, optimcycle,
216     & mythid, xx_sss_dummy )
217     do bj = jtlo,jthi
218     do bi = itlo,ithi
219     do j = jmin,jmax
220     do i = imin,imax
221     cph sss(i,j,bi,bj) = sss(i,j,bi,bj) + tmpfld2d(i,j,bi,bj)
222     salt(i,j,1,bi,bj) = salt(i,j,1,bi,bj)
223     & + tmpfld2d(i,j,bi,bj)
224     enddo
225     enddo
226     enddo
227     enddo
228     #endif
229    
230 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
231     c-- diffkr.
232     il=ilnblnk( xx_diffkr_file )
233 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
234 heimbach 1.3 & xx_diffkr_file(1:il),'.',optimcycle
235 heimbach 1.17 call active_read_xyz_loc( fnamegeneric, tmpfld3d, 1,
236 heimbach 1.3 & doglobalread, ladinit, optimcycle,
237     & mythid, xx_diffkr_dummy )
238     do bj = jtlo,jthi
239     do bi = itlo,ithi
240     do k = 1,nr
241     do j = jmin,jmax
242     do i = imin,imax
243     diffkr(i,j,k,bi,bj) = diffkr(i,j,k,bi,bj) +
244     & tmpfld3d(i,j,k,bi,bj)
245     enddo
246     enddo
247     enddo
248     enddo
249     enddo
250     #endif
251    
252     #ifdef ALLOW_KAPGM_CONTROL
253     c-- kapgm.
254     il=ilnblnk( xx_kapgm_file )
255 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
256 heimbach 1.3 & xx_kapgm_file(1:il),'.',optimcycle
257 heimbach 1.17 call active_read_xyz_loc( fnamegeneric, tmpfld3d, 1,
258 heimbach 1.3 & doglobalread, ladinit, optimcycle,
259     & mythid, xx_kapgm_dummy )
260     do bj = jtlo,jthi
261     do bi = itlo,ithi
262     do k = 1,nr
263     do j = jmin,jmax
264     do i = imin,imax
265     kapgm(i,j,k,bi,bj) = kapgm(i,j,k,bi,bj) +
266     & tmpfld3d(i,j,k,bi,bj)
267     enddo
268     enddo
269     enddo
270     enddo
271     enddo
272     #endif
273    
274 heimbach 1.6 #ifdef ALLOW_EFLUXY0_CONTROL
275     c-- y-component EP-flux field.
276     il=ilnblnk( xx_efluxy_file )
277 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
278 heimbach 1.6 & xx_efluxy_file(1:il),'.',optimcycle
279 heimbach 1.17 call active_read_xyz_loc( fnamegeneric, tmpfld3d, 1,
280 heimbach 1.6 & doglobalread, ladinit, optimcycle,
281     & mythid, xx_efluxy_dummy )
282    
283     do bj = jtlo,jthi
284     do bi = itlo,ithi
285     do k = 1,nr
286     do j = jmin,jmax
287     do i = imin,imax
288     EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
289     & - fac*tmpfld3d(i,j,k,bi,bj)
290     & *maskS(i,j,k,bi,bj)
291     cph EfluxY(i,j,k,bi,bj) = EfluxY(i,j,k,bi,bj)
292     cph & - rSphere*cosFacU(J,bi,bj)
293     cph & *fac*tmpfld3d(i,j,k,bi,bj)
294     enddo
295     enddo
296     enddo
297     enddo
298     enddo
299     #endif
300    
301     #ifdef ALLOW_EFLUXP0_CONTROL
302     c-- p-component EP-flux field.
303     il=ilnblnk( xx_efluxp_file )
304 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
305 heimbach 1.6 & xx_efluxp_file(1:il),'.',optimcycle
306 heimbach 1.17 call active_read_xyz_loc( fnamegeneric, tmpfld3d, 1,
307 heimbach 1.6 & doglobalread, ladinit, optimcycle,
308     & mythid, xx_efluxp_dummy )
309    
310     do bj = jtlo,jthi
311     do bi = itlo,ithi
312     do k = 1,nr
313     do j = jmin,jmax
314     do i = imin,imax
315     EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
316     & + fCori(i,j,bi,bj)
317     & *fac*tmpfld3d(i,j,k,bi,bj)
318     & *hFacV(i,j,k,bi,bj)
319     cph EfluxP(i,j,k,bi,bj) = EfluxP(i,j,k,bi,bj)
320     cph & + fCori(i,j,bi,bj)
321     cph & *rSphere*cosFacU(J,bi,bj)
322     cph & *fac*tmpfld3d(i,j,k,bi,bj)
323     enddo
324     enddo
325     enddo
326     enddo
327     enddo
328     #endif
329    
330 heimbach 1.7 #ifdef ALLOW_BOTTOMDRAG_CONTROL
331     c-- bottom drag
332     il=ilnblnk( xx_bottomdrag_file )
333 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
334 heimbach 1.7 & xx_bottomdrag_file(1:il),'.',optimcycle
335 heimbach 1.17 call active_read_xy_loc ( fnamegeneric, tmpfld2d, 1,
336 heimbach 1.7 & doglobalread, ladinit, optimcycle,
337     & mythid, xx_bottomdrag_dummy )
338     do bj = jtlo,jthi
339     do bi = itlo,ithi
340     do j = jmin,jmax
341     do i = imin,imax
342     bottomdragfld(i,j,bi,bj) = bottomdragfld(i,j,bi,bj)
343     & + tmpfld2d(i,j,bi,bj)
344     enddo
345     enddo
346     enddo
347     enddo
348     #endif
349    
350 heimbach 1.16 #ifdef ALLOW_EDTAUX_CONTROL
351 heimbach 1.15 c-- zonal eddy stress : edtaux
352     il=ilnblnk( xx_edtaux_file )
353 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
354 heimbach 1.15 & xx_edtaux_file(1:il),'.',optimcycle
355 heimbach 1.17 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
356 heimbach 1.15 & doglobalread, ladinit, optimcycle,
357     & mythid, xx_edtaux_dummy )
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     Eddytaux(i,j,k,bi,bj) = Eddytaux(i,j,k,bi,bj) +
364     & tmpfld3d(i,j,k,bi,bj)
365     enddo
366     enddo
367     enddo
368     enddo
369     enddo
370     #endif
371    
372     #ifdef ALLOW_EDTAUY_CONTROL
373     c-- meridional eddy stress : edtauy
374     il=ilnblnk( xx_edtauy_file )
375 heimbach 1.17 write(fnamegeneric(1:80),'(2a,i10.10)')
376 heimbach 1.15 & xx_edtauy_file(1:il),'.',optimcycle
377 heimbach 1.17 call active_read_xyz( fnamegeneric, tmpfld3d, 1,
378 heimbach 1.15 & doglobalread, ladinit, optimcycle,
379     & mythid, xx_edtauy_dummy )
380     do bj = jtlo,jthi
381     do bi = itlo,ithi
382     do k = 1,nr
383     do j = jmin,jmax
384     do i = imin,imax
385     Eddytauy(i,j,k,bi,bj) = Eddytauy(i,j,k,bi,bj) +
386     & tmpfld3d(i,j,k,bi,bj)
387     enddo
388     enddo
389     enddo
390     enddo
391     enddo
392     #endif
393 heimbach 1.1
394 heimbach 1.17 #ifdef ALLOW_UVEL0_CONTROL
395     c-- initial zonal velocity
396     il=ilnblnk( xx_uvel_file )
397     write(fnamegeneric(1:80),'(2a,i10.10)')
398     & xx_uvel_file(1:il),'.',optimcycle
399     call active_read_xyz( fnamegeneric, tmpfld3d, 1,
400     & doglobalread, ladinit, optimcycle,
401     & mythid, xx_uvel_dummy )
402     do bj = jtlo,jthi
403     do bi = itlo,ithi
404     do k = 1,nr
405     do j = jmin,jmax
406     do i = imin,imax
407 heimbach 1.19 #ifdef ALLOW_OPENAD
408     uVel(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) +
409     & fac*xx_uvel(i,j,k,bi,bj)
410     #else
411 heimbach 1.17 uVel(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) +
412 heimbach 1.19 & fac*tmpfld3d(i,j,k,bi,bj)
413     #endif
414 heimbach 1.17 enddo
415     enddo
416     enddo
417     enddo
418     enddo
419     #endif
420    
421     #ifdef ALLOW_VVEL0_CONTROL
422     c-- initial merid. velocity
423     il=ilnblnk( xx_vvel_file )
424     write(fnamegeneric(1:80),'(2a,i10.10)')
425     & xx_vvel_file(1:il),'.',optimcycle
426     call active_read_xyz( fnamegeneric, tmpfld3d, 1,
427     & doglobalread, ladinit, optimcycle,
428     & mythid, xx_vvel_dummy )
429     do bj = jtlo,jthi
430     do bi = itlo,ithi
431     do k = 1,nr
432     do j = jmin,jmax
433     do i = imin,imax
434 heimbach 1.19 #ifdef ALLOW_OPENAD
435 heimbach 1.17 vVel(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) +
436 heimbach 1.19 & fac*xx_vvel(i,j,k,bi,bj)
437     #else
438     vVel(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) +
439     & fac*tmpfld3d(i,j,k,bi,bj)
440     #endif
441 heimbach 1.17 enddo
442     enddo
443     enddo
444     enddo
445     enddo
446     #endif
447    
448     #ifdef ALLOW_ETAN0_CONTROL
449     c-- initial Eta.
450     il=ilnblnk( xx_etan_file )
451     write(fnamegeneric(1:80),'(2a,i10.10)')
452     & xx_etan_file(1:il),'.',optimcycle
453     call active_read_xy_loc ( fnamegeneric, tmpfld2d, 1,
454     & doglobalread, ladinit, optimcycle,
455     & mythid, xx_etan_dummy )
456     do bj = jtlo,jthi
457     do bi = itlo,ithi
458     do j = jmin,jmax
459     do i = imin,imax
460 heimbach 1.19 #ifdef ALLOW_OPENAD
461     etaN(i,j,bi,bj) = etaN(i,j,bi,bj) +
462     & fac*xx_etan(i,j,bi,bj)
463     #else
464     etaN(i,j,bi,bj) = etaN(i,j,bi,bj) +
465     & fac*tmpfld2d(i,j,bi,bj)
466     #endif
467 heimbach 1.17 enddo
468     enddo
469     enddo
470     enddo
471     #endif
472    
473     #ifdef ALLOW_RELAXSST_CONTROL
474     c-- SST relaxation coefficient.
475     il=ilnblnk( xx_relaxsst_file )
476     write(fnamegeneric(1:80),'(2a,i10.10)')
477     & xx_relaxsst_file(1:il),'.',optimcycle
478     call active_read_xy_loc ( fnamegeneric, tmpfld2d, 1,
479     & doglobalread, ladinit, optimcycle,
480     & mythid, xx_relaxsst_dummy )
481     do bj = jtlo,jthi
482     do bi = itlo,ithi
483     do j = jmin,jmax
484     do i = imin,imax
485     lambdaThetaClimRelax(i,j,bi,bj) =
486     & lambdaThetaClimRelax(i,j,bi,bj)
487     & + tmpfld2d(i,j,bi,bj)
488     enddo
489     enddo
490     enddo
491     enddo
492     #endif
493    
494     #ifdef ALLOW_RELAXSSS_CONTROL
495     c-- SSS relaxation coefficient.
496     il=ilnblnk( xx_relaxsss_file )
497     write(fnamegeneric(1:80),'(2a,i10.10)')
498     & xx_relaxsss_file(1:il),'.',optimcycle
499     call active_read_xy_loc ( fnamegeneric, tmpfld2d, 1,
500     & doglobalread, ladinit, optimcycle,
501     & mythid, xx_relaxsss_dummy )
502     do bj = jtlo,jthi
503     do bi = itlo,ithi
504     do j = jmin,jmax
505     do i = imin,imax
506     lambdaSaltClimRelax(i,j,bi,bj) =
507     & lambdaSaltClimRelax(i,j,bi,bj)
508     & + tmpfld2d(i,j,bi,bj)
509     enddo
510     enddo
511     enddo
512     enddo
513     #endif
514    
515 heimbach 1.1 c-- Update the tile edges.
516    
517 heimbach 1.14 #if (defined (ALLOW_THETA0_CONTROL) || defined (ALLOW_SST0_CONTROL))
518 heimbach 1.1 _EXCH_XYZ_R8( theta, mythid )
519     #endif
520 heimbach 1.14 #if (defined (ALLOW_SALT0_CONTROL) || defined (ALLOW_SSS0_CONTROL))
521 heimbach 1.1 _EXCH_XYZ_R8( salt, mythid )
522 heimbach 1.2 #endif
523     #ifdef ALLOW_TR10_CONTROL
524 heimbach 1.13 #ifdef ALLOW_PTRACERS
525     _EXCH_XYZ_R8(pTracer(1-Olx,1-Oly,1,1,1,1),myThid)
526     #endif
527 heimbach 1.1 #endif
528 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
529     _EXCH_XYZ_R8( diffkr, mythid)
530     #endif
531     #ifdef ALLOW_KAPGM_CONTROL
532     _EXCH_XYZ_R8( kapgm, mythid)
533 heimbach 1.6 #endif
534     #ifdef ALLOW_EFLUXY0_CONTROL
535     _EXCH_XYZ_R8( EfluxY, mythid )
536     #endif
537     #ifdef ALLOW_EFLUXP0_CONTROL
538     _EXCH_XYZ_R8( EfluxP, mythid )
539 heimbach 1.7 #endif
540     #ifdef ALLOW_BOTTOMDRAG_CONTROL
541     _EXCH_XY_R8( bottomdragfld, mythid )
542 heimbach 1.3 #endif
543    
544 heimbach 1.15 #if (defined (ALLOW_EDTAUX_CONTROL) && defined (ALLOW_EDTAUY_CONTROL))
545     CALL EXCH_UV_XYZ_RS(Eddytaux,Eddytauy,.TRUE.,myThid)
546     #elif (defined (ALLOW_EDTAUX_CONTROL) || defined (ALLOW_EDTAUY_CONTROL))
547     STOP 'ctrl_map_forcing: need BOTH ALLOW_EDTAU[X,Y]_CONTROL'
548     #endif
549 heimbach 1.1
550 heimbach 1.17 #ifdef ALLOW_UVEL0_CONTROL
551     _EXCH_XYZ_R8( uVel, mythid)
552     #endif
553    
554     #ifdef ALLOW_VVEL0_CONTROL
555     _EXCH_XYZ_R8( vVel, mythid)
556     #endif
557    
558     #ifdef ALLOW_ETAN0_CONTROL
559     _EXCH_XY_R8( etaN, mythid )
560     #endif
561    
562     #ifdef ALLOW_RELAXSST_CONTROL
563     _EXCH_XY_R4( lambdaThetaClimRelax, mythid )
564     #endif
565    
566     #ifdef ALLOW_RELAXSSS_CONTROL
567     _EXCH_XY_R4( lambdaThetaClimRelax, mythid )
568     #endif
569    
570 heimbach 1.1 return
571     end
572    

  ViewVC Help
Powered by ViewVC 1.1.22