/[MITgcm]/MITgcm/pkg/ctrl/ctrl_map_ini.F
ViewVC logotype

Diff of /MITgcm/pkg/ctrl/ctrl_map_ini.F

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

revision 1.10 by heimbach, Fri Jun 27 01:54:20 2003 UTC revision 1.31 by utke, Tue Nov 18 16:47:15 2008 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2    C $Name$
3    
4  #include "CTRL_CPPOPTIONS.h"  #include "CTRL_CPPOPTIONS.h"
 #ifdef ALLOW_PTRACERS  
 # include "PTRACERS_OPTIONS.h"  
 #endif  
5    
6  CBOP  CBOP
7  C     !ROUTINE: ctrl_map_ini  C     !ROUTINE: ctrl_map_ini
# Line 13  C     !INTERFACE: Line 11  C     !INTERFACE:
11  C     !DESCRIPTION: \bv  C     !DESCRIPTION: \bv
12  c     *=================================================================  c     *=================================================================
13  c     | SUBROUTINE ctrl_map_ini  c     | SUBROUTINE ctrl_map_ini
14  c     | Add the temperature, salinity, and diffusivity parts of the  c     | Add the temperature, salinity, and diffusivity parts of the
15  c     | control vector to the model state and update the tile halos.  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".  c     | The control vector is defined in the header file "ctrl.h".
17  c     *=================================================================  c     *=================================================================
18  C     \ev  C     \ev
# Line 26  c     == global variables == Line 24  c     == global variables ==
24  #include "SIZE.h"  #include "SIZE.h"
25  #include "EEPARAMS.h"  #include "EEPARAMS.h"
26  #include "PARAMS.h"  #include "PARAMS.h"
 #include "DYNVARS.h"  
27  #include "GRID.h"  #include "GRID.h"
28  #ifdef ALLOW_PASSIVE_TRACER  #include "DYNVARS.h"
29  # include "TR1.h"  #include "FFIELDS.h"
 #endif  
 #ifdef ALLOW_PTRACERS  
 # include "PTRACERS.h"  
 #endif  
   
30  #include "ctrl.h"  #include "ctrl.h"
31  #include "ctrl_dummy.h"  #include "ctrl_dummy.h"
32  #include "optim.h"  #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    
42  C     !INPUT/OUTPUT PARAMETERS:  C     !INPUT/OUTPUT PARAMETERS:
43  c     == routine arguments ==  c     == routine arguments ==
# Line 58  c     == local variables == Line 58  c     == local variables ==
58        logical doglobalread        logical doglobalread
59        logical ladinit        logical ladinit
60    
61        character*( 80)   fnametheta        character*( 80)   fnamegeneric
       character*( 80)   fnamesalt  
       character*( 80)   fnametr1  
       character*( 80)   fnamediffkr  
       character*( 80)   fnamekapgm  
       character*( 80)   fnameefluxy  
       character*( 80)   fnameefluxp  
       character*( 80)   fnamebottomdrag  
62    
63        _RL     fac        _RL     fac
64          _RL tmptest
65    
66  c     == external ==  c     == external ==
67        integer  ilnblnk        integer  ilnblnk
# Line 99  CEOP Line 93  CEOP
93  #ifdef ALLOW_THETA0_CONTROL  #ifdef ALLOW_THETA0_CONTROL
94  c--   Temperature field.  c--   Temperature field.
95        il=ilnblnk( xx_theta_file )        il=ilnblnk( xx_theta_file )
96        write(fnametheta(1:80),'(2a,i10.10)')        write(fnamegeneric(1:80),'(2a,i10.10)')
97       &     xx_theta_file(1:il),'.',optimcycle       &     xx_theta_file(1:il),'.',optimcycle
98        call active_read_xyz( fnametheta, tmpfld3d, 1,        call active_read_xyz( fnamegeneric, tmpfld3d, 1,
99       &                      doglobalread, ladinit, optimcycle,       &                      doglobalread, ladinit, optimcycle,
100       &                      mythid, xx_theta_dummy )       &                      mythid, xx_theta_dummy )
101    
# Line 110  c--   Temperature field. Line 104  c--   Temperature field.
104            do k = 1,nr            do k = 1,nr
105              do j = jmin,jmax              do j = jmin,jmax
106                do i = imin,imax                do i = imin,imax
107    #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    #ifdef ALLOW_AUTODIFF_OPENAD
114                  theta(i,j,k,bi,bj) = theta(i,j,k,bi,bj) +                  theta(i,j,k,bi,bj) = theta(i,j,k,bi,bj) +
115         &                               fac*xx_theta(i,j,k,bi,bj) +
116       &                               fac*tmpfld3d(i,j,k,bi,bj)       &                               fac*tmpfld3d(i,j,k,bi,bj)
117                  if(theta(i,j,k,bi,bj).lt.-2.0)  #else
118       &               theta(i,j,k,bi,bj)= -2.0                    theta(i,j,k,bi,bj) = theta(i,j,k,bi,bj) +
119         &                               fac*tmpfld3d(i,j,k,bi,bj)
120    #endif
121    #ifndef DISABLE_CTRL_THETA_LIMIT
122                    if(theta(i,j,k,bi,bj).lt.-2.0)
123         &               theta(i,j,k,bi,bj)= -2.0
124    #endif
125                enddo                enddo
126              enddo              enddo
127            enddo            enddo
128         enddo         enddo
129        enddo        enddo
130    
131  #endif  #endif
132    
133  #ifdef ALLOW_SALT0_CONTROL  #ifdef ALLOW_SALT0_CONTROL
134  c--   Temperature field.  c--   Temperature field.
135        il=ilnblnk( xx_salt_file )        il=ilnblnk( xx_salt_file )
136        write(fnamesalt(1:80),'(2a,i10.10)')        write(fnamegeneric(1:80),'(2a,i10.10)')
137       &     xx_salt_file(1:il),'.',optimcycle       &     xx_salt_file(1:il),'.',optimcycle
138        call active_read_xyz( fnamesalt, tmpfld3d, 1,        call active_read_xyz( fnamegeneric, tmpfld3d, 1,
139       &                      doglobalread, ladinit, optimcycle,       &                      doglobalread, ladinit, optimcycle,
140       &                      mythid, xx_salt_dummy )       &                      mythid, xx_salt_dummy )
141    
# Line 135  c--   Temperature field. Line 144  c--   Temperature field.
144            do k = 1,nr            do k = 1,nr
145              do j = jmin,jmax              do j = jmin,jmax
146                do i = imin,imax                do i = imin,imax
147    #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    #ifdef ALLOW_AUTODIFF_OPENAD
154                  salt(i,j,k,bi,bj) = salt(i,j,k,bi,bj) +                  salt(i,j,k,bi,bj) = salt(i,j,k,bi,bj) +
155         &                               fac*xx_salt(i,j,k,bi,bj) +
156       &                               fac*tmpfld3d(i,j,k,bi,bj)       &                               fac*tmpfld3d(i,j,k,bi,bj)
157    #else
158                    salt(i,j,k,bi,bj) = salt(i,j,k,bi,bj) +
159         &                               fac*tmpfld3d(i,j,k,bi,bj)
160    #endif
161    
162                enddo                enddo
163              enddo              enddo
164            enddo            enddo
165         enddo         enddo
166        enddo        enddO
167  #endif  #endif
168    
169  #ifdef ALLOW_TR10_CONTROL  #ifdef ALLOW_TR10_CONTROL
170    #ifdef ALLOW_PTRACERS
171  c--   Temperature field.  c--   Temperature field.
172        il=ilnblnk( xx_tr1_file )        il=ilnblnk( xx_tr1_file )
173        write(fnametr1(1:80),'(2a,i10.10)')        write(fnamegeneric(1:80),'(2a,i10.10)')
174       &     xx_tr1_file(1:il),'.',optimcycle       &     xx_tr1_file(1:il),'.',optimcycle
175        call active_read_xyz( fnametr1, tmpfld3d, 1,        call active_read_xyz( fnamegeneric, tmpfld3d, 1,
176       &                      doglobalread, ladinit, optimcycle,       &                      doglobalread, ladinit, optimcycle,
177       &                      mythid, xx_tr1_dummy )       &                      mythid, xx_tr1_dummy )
178    
# Line 158  c--   Temperature field. Line 181  c--   Temperature field.
181            do k = 1,nr            do k = 1,nr
182              do j = jmin,jmax              do j = jmin,jmax
183                do i = imin,imax                do i = imin,imax
 #if (defined (ALLOW_PASSIVE_TRACER))  
                 tr1(i,j,k,bi,bj) = tr1(i,j,k,bi,bj) +  
      &                               fac*tmpfld3d(i,j,k,bi,bj)  
 #elif (defined (ALLOW_PTRACERS))  
                 IF ( NUMBER_OF_PTRACERS .GT. 1 ) STOP  
      & 'ALLOW_TR10_CONTROL with ALLOW_PTRACERS implemented for 1 tracer'  
184                  ptracer(i,j,k,bi,bj,1) = ptracer(i,j,k,bi,bj,1) +                  ptracer(i,j,k,bi,bj,1) = ptracer(i,j,k,bi,bj,1) +
185       &                               fac*tmpfld3d(i,j,k,bi,bj)       &                               fac*tmpfld3d(i,j,k,bi,bj)
 #endif  
186                enddo                enddo
187              enddo              enddo
188            enddo            enddo
189         enddo         enddo
190        enddo        enddo
191  #endif  #endif
192    #endif
193    
194    #ifdef ALLOW_SST0_CONTROL
195    c--   sst0.
196          il=ilnblnk( xx_sst_file )
197          write(fnamegeneric(1:80),'(2a,i10.10)')
198         &     xx_sst_file(1:il),'.',optimcycle
199          call active_read_xy ( fnamegeneric, tmpfld2d, 1,
200         &                      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                  theta(i,j,1,bi,bj) = theta(i,j,1,bi,bj)
208         &                             + 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          write(fnamegeneric(1:80),'(2a,i10.10)')
219         &     xx_sss_file(1:il),'.',optimcycle
220          call active_read_xy ( fnamegeneric, tmpfld2d, 1,
221         &                      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    #ifdef ALLOW_AUTODIFF
237  #ifdef ALLOW_DIFFKR_CONTROL  #ifdef ALLOW_DIFFKR_CONTROL
238  c--   diffkr.  c--   diffkr.
239        il=ilnblnk( xx_diffkr_file )        il=ilnblnk( xx_diffkr_file )
240        write(fnamediffkr(1:80),'(2a,i10.10)')        write(fnamegeneric(1:80),'(2a,i10.10)')
241       &     xx_diffkr_file(1:il),'.',optimcycle       &     xx_diffkr_file(1:il),'.',optimcycle
242        call active_read_xyz( fnamediffkr, tmpfld3d, 1,        call active_read_xyz( fnamegeneric, tmpfld3d, 1,
243       &                      doglobalread, ladinit, optimcycle,       &                      doglobalread, ladinit, optimcycle,
244       &                      mythid, xx_diffkr_dummy )       &                      mythid, xx_diffkr_dummy )
245        do bj = jtlo,jthi        do bj = jtlo,jthi
# Line 187  c--   diffkr. Line 247  c--   diffkr.
247            do k = 1,nr            do k = 1,nr
248              do j = jmin,jmax              do j = jmin,jmax
249                do i = imin,imax                do i = imin,imax
250    #ifdef ALLOW_AUTODIFF_OPENAD
251                  diffkr(i,j,k,bi,bj) = diffkr(i,j,k,bi,bj) +                  diffkr(i,j,k,bi,bj) = diffkr(i,j,k,bi,bj) +
252         &                                xx_diffkr(i,j,k,bi,bj) +
253       &                                tmpfld3d(i,j,k,bi,bj)       &                                tmpfld3d(i,j,k,bi,bj)
254    #else
255                    diffkr(i,j,k,bi,bj) = diffkr(i,j,k,bi,bj) +
256         &                                tmpfld3d(i,j,k,bi,bj)
257    #endif
258                enddo                enddo
259              enddo              enddo
260            enddo            enddo
261         enddo         enddo
262        enddo        enddo
263  #endif  #endif
264    #endif
265    
266    #ifdef ALLOW_AUTODIFF
267  #ifdef ALLOW_KAPGM_CONTROL  #ifdef ALLOW_KAPGM_CONTROL
268  c--   kapgm.  c--   kapgm.
269        il=ilnblnk( xx_kapgm_file )        il=ilnblnk( xx_kapgm_file )
270        write(fnamekapgm(1:80),'(2a,i10.10)')        write(fnamegeneric(1:80),'(2a,i10.10)')
271       &     xx_kapgm_file(1:il),'.',optimcycle       &     xx_kapgm_file(1:il),'.',optimcycle
272        call active_read_xyz( fnamekapgm, tmpfld3d, 1,        call active_read_xyz( fnamegeneric, tmpfld3d, 1,
273       &                      doglobalread, ladinit, optimcycle,       &                      doglobalread, ladinit, optimcycle,
274       &                      mythid, xx_kapgm_dummy )       &                      mythid, xx_kapgm_dummy )
275        do bj = jtlo,jthi        do bj = jtlo,jthi
# Line 209  c--   kapgm. Line 277  c--   kapgm.
277            do k = 1,nr            do k = 1,nr
278              do j = jmin,jmax              do j = jmin,jmax
279                do i = imin,imax                do i = imin,imax
280    #ifdef ALLOW_AUTODIFF_OPENAD
281                    kapgm(i,j,k,bi,bj) = kapgm(i,j,k,bi,bj) +
282         &                               xx_kapgm(i,j,k,bi,bj) +
283         &                               tmpfld3d(i,j,k,bi,bj)
284    #else
285                  kapgm(i,j,k,bi,bj) = kapgm(i,j,k,bi,bj) +                  kapgm(i,j,k,bi,bj) = kapgm(i,j,k,bi,bj) +
286       &                               tmpfld3d(i,j,k,bi,bj)       &                               tmpfld3d(i,j,k,bi,bj)
287    #endif
288                  enddo
289                enddo
290              enddo
291           enddo
292          enddo
293    #endif
294    #endif
295    
296    #ifdef ALLOW_AUTODIFF
297    #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                enddo
313              enddo              enddo
314            enddo            enddo
315         enddo         enddo
316        enddo        enddo
317  #endif  #endif
318    #endif
319    
320  #ifdef ALLOW_EFLUXY0_CONTROL  #ifdef ALLOW_EFLUXY0_CONTROL
321  c--   y-component EP-flux field.  c--   y-component EP-flux field.
322        il=ilnblnk( xx_efluxy_file )        il=ilnblnk( xx_efluxy_file )
323        write(fnameefluxy(1:80),'(2a,i10.10)')        write(fnamegeneric(1:80),'(2a,i10.10)')
324       &     xx_efluxy_file(1:il),'.',optimcycle       &     xx_efluxy_file(1:il),'.',optimcycle
325        call active_read_xyz( fnameefluxy, tmpfld3d, 1,        call active_read_xyz( fnamegeneric, tmpfld3d, 1,
326       &                      doglobalread, ladinit, optimcycle,       &                      doglobalread, ladinit, optimcycle,
327       &                      mythid, xx_efluxy_dummy )       &                      mythid, xx_efluxy_dummy )
328    
# Line 248  cph     & Line 347  cph     &
347  #ifdef ALLOW_EFLUXP0_CONTROL  #ifdef ALLOW_EFLUXP0_CONTROL
348  c--   p-component EP-flux field.  c--   p-component EP-flux field.
349        il=ilnblnk( xx_efluxp_file )        il=ilnblnk( xx_efluxp_file )
350        write(fnameefluxp(1:80),'(2a,i10.10)')        write(fnamegeneric(1:80),'(2a,i10.10)')
351       &     xx_efluxp_file(1:il),'.',optimcycle       &     xx_efluxp_file(1:il),'.',optimcycle
352        call active_read_xyz( fnameefluxp, tmpfld3d, 1,        call active_read_xyz( fnamegeneric, tmpfld3d, 1,
353       &                      doglobalread, ladinit, optimcycle,       &                      doglobalread, ladinit, optimcycle,
354       &                      mythid, xx_efluxp_dummy )       &                      mythid, xx_efluxp_dummy )
355    
# Line 277  cph     & Line 376  cph     &
376  #ifdef ALLOW_BOTTOMDRAG_CONTROL  #ifdef ALLOW_BOTTOMDRAG_CONTROL
377  c--   bottom drag  c--   bottom drag
378        il=ilnblnk( xx_bottomdrag_file )        il=ilnblnk( xx_bottomdrag_file )
379        write(fnamebottomdrag(1:80),'(2a,i10.10)')        write(fnamegeneric(1:80),'(2a,i10.10)')
380       &     xx_bottomdrag_file(1:il),'.',optimcycle       &     xx_bottomdrag_file(1:il),'.',optimcycle
381        call active_read_xy ( fnamebottomdrag, tmpfld2d, 1,        call active_read_xy ( fnamegeneric, tmpfld2d, 1,
382       &                      doglobalread, ladinit, optimcycle,       &                      doglobalread, ladinit, optimcycle,
383       &                      mythid, xx_bottomdrag_dummy )       &                      mythid, xx_bottomdrag_dummy )
384        do bj = jtlo,jthi        do bj = jtlo,jthi
385          do bi = itlo,ithi          do bi = itlo,ithi
386            do j = jmin,jmax            do j = jmin,jmax
387              do i = imin,imax              do i = imin,imax
388                bottomdragfld(i,j,bi,bj) = bottomdragfld(i,j,bi,bj)                bottomdragfld(i,j,bi,bj) = bottomdragfld(i,j,bi,bj)
389       &                                   + tmpfld2d(i,j,bi,bj)       &                                   + tmpfld2d(i,j,bi,bj)
390              enddo              enddo
391            enddo            enddo
# Line 294  c--   bottom drag Line 393  c--   bottom drag
393        enddo        enddo
394  #endif  #endif
395    
396    #ifdef ALLOW_EDDYPSI_CONTROL
397    c-- zonal eddy streamfunction : eddyPsiX
398          il=ilnblnk( xx_edtaux_file )
399          write(fnamegeneric(1:80),'(2a,i10.10)')
400         &     xx_edtaux_file(1:il),'.',optimcycle
401          call active_read_xyz( fnamegeneric, tmpfld3d, 1,
402         &                      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                    eddyPsiX(i,j,k,bi,bj) = eddyPsiX(i,j,k,bi,bj) +
410         &            tmpfld3d(i,j,k,bi,bj)
411                  enddo
412                enddo
413              enddo
414           enddo
415          enddo
416    c-- meridional eddy streamfunction : eddyPsiY
417          il=ilnblnk( xx_edtauy_file )
418          write(fnamegeneric(1:80),'(2a,i10.10)')
419         &     xx_edtauy_file(1:il),'.',optimcycle
420          call active_read_xyz( fnamegeneric, tmpfld3d, 1,
421         &                      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                    eddyPsiY(i,j,k,bi,bj) = eddyPsiY(i,j,k,bi,bj) +
429         &            tmpfld3d(i,j,k,bi,bj)
430                  enddo
431                enddo
432              enddo
433           enddo
434          enddo
435    #endif
436    
437    #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    #ifdef ALLOW_AUTODIFF_OPENAD
451                    uVel(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) +
452         &                                  fac*xx_uvel(i,j,k,bi,bj)
453    #else
454                    uVel(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) +
455         &                                  fac*tmpfld3d(i,j,k,bi,bj)
456    #endif
457                  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    #ifdef ALLOW_AUTODIFF_OPENAD
478                    vVel(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) +
479         &                                  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                  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          call active_read_xy ( fnamegeneric, tmpfld2d, 1,
497         &                      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    #ifdef ALLOW_AUTODIFF_OPENAD
504                  etaN(i,j,bi,bj) = etaN(i,j,bi,bj) +
505         &                              fac*xx_etan(i,j,bi,bj)
506    #else
507                  etaN(i,j,bi,bj) = etaN(i,j,bi,bj) +
508         &                              fac*tmpfld2d(i,j,bi,bj)
509    #endif
510                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          call active_read_xy ( fnamegeneric, tmpfld2d, 1,
522         &                      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                  lambdaThetaClimRelax(i,j,bi,bj) =
529         &              lambdaThetaClimRelax(i,j,bi,bj)
530         &              + 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          call active_read_xy ( fnamegeneric, tmpfld2d, 1,
543         &                      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                  lambdaSaltClimRelax(i,j,bi,bj) =
550         &              lambdaSaltClimRelax(i,j,bi,bj)
551         &              + tmpfld2d(i,j,bi,bj)
552                enddo
553              enddo
554            enddo
555          enddo
556    #endif
557    
558    #ifdef ALLOW_SEAICE
559          call seaice_ctrl_map_ini( mythid )
560    #endif
561    
562  c--   Update the tile edges.  c--   Update the tile edges.
563    
564  #ifdef ALLOW_THETA0_CONTROL  #if (defined (ALLOW_THETA0_CONTROL) || defined (ALLOW_SST0_CONTROL))
565        _EXCH_XYZ_R8( theta, mythid )        _EXCH_XYZ_R8( theta, mythid )
566  #endif  #endif
567  #ifdef ALLOW_SALT0_CONTROL  #if (defined (ALLOW_SALT0_CONTROL) || defined (ALLOW_SSS0_CONTROL))
568        _EXCH_XYZ_R8(  salt, mythid )        _EXCH_XYZ_R8(  salt, mythid )
569  #endif  #endif
570  #ifdef ALLOW_TR10_CONTROL  #ifdef ALLOW_TR10_CONTROL
571  # if (defined (ALLOW_PASSIVE_TRACER))  #ifdef ALLOW_PTRACERS
       _EXCH_XYZ_R8(     tr1, mythid )  
 # elif (defined (ALLOW_PTRACERS))  
572        _EXCH_XYZ_R8(pTracer(1-Olx,1-Oly,1,1,1,1),myThid)        _EXCH_XYZ_R8(pTracer(1-Olx,1-Oly,1,1,1,1),myThid)
 # endif  
573  #endif  #endif
 #ifdef ALLOW_DIFFKR_CONTROL  
       _EXCH_XYZ_R8( diffkr, mythid)  
574  #endif  #endif
575  #ifdef ALLOW_KAPGM_CONTROL  
576    #ifdef ALLOW_AUTODIFF
577    # ifdef ALLOW_DIFFKR_CONTROL
578          _EXCH_XYZ_R8( diffkr, mythid)
579    # endif
580    # ifdef ALLOW_KAPGM_CONTROL
581        _EXCH_XYZ_R8( kapgm, mythid)        _EXCH_XYZ_R8( kapgm, mythid)
582    # endif
583    # ifdef ALLOW_KAPREDI_CONTROL
584          _EXCH_XYZ_R8( kapredi, mythid)
585    # endif
586  #endif  #endif
587    
588  #ifdef ALLOW_EFLUXY0_CONTROL  #ifdef ALLOW_EFLUXY0_CONTROL
589        _EXCH_XYZ_R8( EfluxY, mythid )        _EXCH_XYZ_R8( EfluxY, mythid )
590  #endif  #endif
# Line 326  c--   Update the tile edges. Line 595  c--   Update the tile edges.
595        _EXCH_XY_R8( bottomdragfld, mythid )        _EXCH_XY_R8( bottomdragfld, mythid )
596  #endif  #endif
597    
598    #ifdef ALLOW_EDDYPSI_CONTROL
599           CALL EXCH_UV_XYZ_RS(eddyPsiX,eddyPsiY,.TRUE.,myThid)
600    #endif
601    
602    #ifdef ALLOW_UVEL0_CONTROL
603          _EXCH_XYZ_R8( uVel, mythid)
604    #endif
605    
606    #ifdef ALLOW_VVEL0_CONTROL
607          _EXCH_XYZ_R8( vVel, mythid)
608    #endif
609    
610    #ifdef ALLOW_ETAN0_CONTROL
611          _EXCH_XY_R8( etaN, mythid )
612    #endif
613    
614    #ifdef ALLOW_RELAXSST_CONTROL
615          _EXCH_XY_R4( lambdaThetaClimRelax, mythid )
616    #endif
617    
618    #ifdef ALLOW_RELAXSSS_CONTROL
619          _EXCH_XY_R4( lambdaThetaClimRelax, mythid )
620    #endif
621    
622        return        return
623        end        end

Legend:
Removed from v.1.10  
changed lines
  Added in v.1.31

  ViewVC Help
Powered by ViewVC 1.1.22