New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
traadv.F90 in trunk/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: trunk/NEMOGCM/NEMO/OPA_SRC/TRA/traadv.F90 @ 4624

Last change on this file since 4624 was 4624, checked in by acc, 10 years ago

#1305. Fix slow start-up problems on some systems by introducing and using lwm logical to restrict output of merged namelists to the first (or only) processor. lwm is true only on the first processor regardless of ln_ctl. Small changes to all flavours of nemogcm.F90 are also required to write namctl and namcfg after the call to mynode which now opens output.namelist.dyn and writes nammpp.

  • Property svn:keywords set to Id
File size: 13.1 KB
Line 
1MODULE traadv
2   !!==============================================================================
3   !!                       ***  MODULE  traadv  ***
4   !! Ocean active tracers:  advection trend
5   !!==============================================================================
6   !! History :  2.0  !  2005-11  (G. Madec)  Original code
7   !!            3.3  !  2010-09  (C. Ethe, G. Madec)  merge TRC-TRA + switch from velocity to transport
8   !!            4.0  !  2011-06  (G. Madec)  Addition of Mixed Layer Eddy parameterisation
9   !!----------------------------------------------------------------------
10
11   !!----------------------------------------------------------------------
12   !!   tra_adv      : compute ocean tracer advection trend
13   !!   tra_adv_ctl  : control the different options of advection scheme
14   !!----------------------------------------------------------------------
15   USE oce             ! ocean dynamics and active tracers
16   USE dom_oce         ! ocean space and time domain
17   USE domvvl          ! variable vertical scale factors
18   USE traadv_cen2     ! 2nd order centered scheme (tra_adv_cen2   routine)
19   USE traadv_tvd      ! TVD      scheme           (tra_adv_tvd    routine)
20   USE traadv_muscl    ! MUSCL    scheme           (tra_adv_muscl  routine)
21   USE traadv_muscl2   ! MUSCL2   scheme           (tra_adv_muscl2 routine)
22   USE traadv_ubs      ! UBS      scheme           (tra_adv_ubs    routine)
23   USE traadv_qck      ! QUICKEST scheme           (tra_adv_qck    routine)
24   USE traadv_eiv      ! eddy induced velocity     (tra_adv_eiv    routine)
25   USE traadv_mle      ! ML eddy induced velocity  (tra_adv_mle    routine)
26   USE cla             ! cross land advection      (cla_traadv     routine)
27   USE ldftra_oce      ! lateral diffusion coefficient on tracers
28   USE in_out_manager  ! I/O manager
29   USE iom             ! I/O module
30   USE prtctl          ! Print control
31   USE lib_mpp         ! MPP library
32   USE wrk_nemo        ! Memory Allocation
33   USE timing          ! Timing
34
35
36   IMPLICIT NONE
37   PRIVATE
38
39   PUBLIC   tra_adv        ! routine called by step module
40   PUBLIC   tra_adv_init   ! routine called by opa module
41
42   !                              !!* Namelist namtra_adv *
43   LOGICAL ::   ln_traadv_cen2     ! 2nd order centered scheme flag
44   LOGICAL ::   ln_traadv_tvd      ! TVD scheme flag
45   LOGICAL ::   ln_traadv_muscl    ! MUSCL scheme flag
46   LOGICAL ::   ln_traadv_muscl2   ! MUSCL2 scheme flag
47   LOGICAL ::   ln_traadv_ubs      ! UBS scheme flag
48   LOGICAL ::   ln_traadv_qck      ! QUICKEST scheme flag
49   LOGICAL ::   ln_traadv_msc_ups  ! use upstream scheme within muscl
50
51
52   INTEGER ::   nadv   ! choice of the type of advection scheme
53
54   !! * Substitutions
55#  include "domzgr_substitute.h90"
56#  include "vectopt_loop_substitute.h90"
57   !!----------------------------------------------------------------------
58   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
59   !! $Id$
60   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
61   !!----------------------------------------------------------------------
62CONTAINS
63
64   SUBROUTINE tra_adv( kt )
65      !!----------------------------------------------------------------------
66      !!                  ***  ROUTINE tra_adv  ***
67      !!
68      !! ** Purpose :   compute the ocean tracer advection trend.
69      !!
70      !! ** Method  : - Update (ua,va) with the advection term following nadv
71      !!----------------------------------------------------------------------
72      !
73      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
74      !
75      INTEGER ::   jk   ! dummy loop index
76      REAL(wp), POINTER, DIMENSION(:,:,:) :: zun, zvn, zwn
77      !!----------------------------------------------------------------------
78      !
79      IF( nn_timing == 1 )  CALL timing_start('tra_adv')
80      !
81      CALL wrk_alloc( jpi, jpj, jpk, zun, zvn, zwn )
82      !                                          ! set time step
83      IF( neuler == 0 .AND. kt == nit000 ) THEN     ! at nit000
84         r2dtra(:) =  rdttra(:)                          ! = rdtra (restarting with Euler time stepping)
85      ELSEIF( kt <= nit000 + 1) THEN                ! at nit000 or nit000+1
86         r2dtra(:) = 2._wp * rdttra(:)                   ! = 2 rdttra (leapfrog)
87      ENDIF
88      !
89      IF( nn_cla == 1 .AND. cp_cfg == 'orca' .AND. jp_cfg == 2 )   CALL cla_traadv( kt )       !==  Cross Land Advection  ==! (hor. advection)
90      !
91      !                                               !==  effective transport  ==!
92      DO jk = 1, jpkm1
93         zun(:,:,jk) = e2u(:,:) * fse3u(:,:,jk) * un(:,:,jk)                  ! eulerian transport only
94         zvn(:,:,jk) = e1v(:,:) * fse3v(:,:,jk) * vn(:,:,jk)
95         zwn(:,:,jk) = e1t(:,:) * e2t(:,:)      * wn(:,:,jk)
96      END DO
97      !
98      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
99         zun(:,:,:) = zun(:,:,:) + un_td(:,:,:)
100         zvn(:,:,:) = zvn(:,:,:) + vn_td(:,:,:)
101      ENDIF
102      !
103      zun(:,:,jpk) = 0._wp                                                     ! no transport trough the bottom
104      zvn(:,:,jpk) = 0._wp                                                     ! no transport trough the bottom
105      zwn(:,:,jpk) = 0._wp                                                     ! no transport trough the bottom
106      !
107      IF( lk_traldf_eiv .AND. .NOT. ln_traldf_grif )   &
108         &              CALL tra_adv_eiv( kt, nit000, zun, zvn, zwn, 'TRA' )    ! add the eiv transport (if necessary)
109      !
110      IF( ln_mle    )   CALL tra_adv_mle( kt, nit000, zun, zvn, zwn, 'TRA' )    ! add the mle transport (if necessary)
111      CALL iom_put( "uocetr_eff", zun )                                         ! output effective transport     
112      CALL iom_put( "vocetr_eff", zvn )
113      CALL iom_put( "wocetr_eff", zwn )
114
115      SELECT CASE ( nadv )                            !==  compute advection trend and add it to general trend  ==!
116      CASE ( 1 )   ;    CALL tra_adv_cen2  ( kt, nit000, 'TRA',         zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  2nd order centered
117      CASE ( 2 )   ;    CALL tra_adv_tvd   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  TVD
118      CASE ( 3 )   ;    CALL tra_adv_muscl ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb,      tsa, jpts, ln_traadv_msc_ups )   !  MUSCL
119      CASE ( 4 )   ;    CALL tra_adv_muscl2( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  MUSCL2
120      CASE ( 5 )   ;    CALL tra_adv_ubs   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  UBS
121      CASE ( 6 )   ;    CALL tra_adv_qck   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  QUICKEST
122      !
123      CASE (-1 )                                      !==  esopa: test all possibility with control print  ==!
124         CALL tra_adv_cen2  ( kt, nit000, 'TRA',         zun, zvn, zwn, tsb, tsn, tsa, jpts )         
125         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv0 - Ta: ', mask1=tmask,               &
126            &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
127         CALL tra_adv_tvd   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )         
128         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv1 - Ta: ', mask1=tmask,               &
129            &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
130         CALL tra_adv_muscl ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb,      tsa, jpts, ln_traadv_msc_ups )         
131         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv3 - Ta: ', mask1=tmask,               &
132            &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
133         CALL tra_adv_muscl2( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )         
134         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv4 - Ta: ', mask1=tmask,               &
135            &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
136         CALL tra_adv_ubs   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )         
137         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv5 - Ta: ', mask1=tmask,               &
138            &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
139         CALL tra_adv_qck   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )         
140         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv6 - Ta: ', mask1=tmask,               &
141            &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
142      END SELECT
143      !
144      !                                              ! print mean trends (used for debugging)
145      IF(ln_ctl)   CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv  - Ta: ', mask1=tmask,               &
146         &                       tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
147      !
148      IF( nn_timing == 1 )  CALL timing_stop( 'tra_adv' )
149      !
150      CALL wrk_dealloc( jpi, jpj, jpk, zun, zvn, zwn )
151      !                                         
152   END SUBROUTINE tra_adv
153
154
155   SUBROUTINE tra_adv_init
156      !!---------------------------------------------------------------------
157      !!                  ***  ROUTINE tra_adv_init  ***
158      !!               
159      !! ** Purpose :   Control the consistency between namelist options for
160      !!              tracer advection schemes and set nadv
161      !!----------------------------------------------------------------------
162      INTEGER ::   ioptio
163      INTEGER ::   ios                 ! Local integer output status for namelist read
164      !!
165      NAMELIST/namtra_adv/ ln_traadv_cen2 , ln_traadv_tvd,     &
166         &                 ln_traadv_muscl, ln_traadv_muscl2,  &
167         &                 ln_traadv_ubs  , ln_traadv_qck,     &
168         &                 ln_traadv_msc_ups
169      !!----------------------------------------------------------------------
170
171      REWIND( numnam_ref )              ! Namelist namtra_adv in reference namelist : Tracer advection scheme
172      READ  ( numnam_ref, namtra_adv, IOSTAT = ios, ERR = 901)
173901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_adv in reference namelist', lwp )
174
175      REWIND( numnam_cfg )              ! Namelist namtra_adv in configuration namelist : Tracer advection scheme
176      READ  ( numnam_cfg, namtra_adv, IOSTAT = ios, ERR = 902 )
177902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_adv in configuration namelist', lwp )
178      IF(lwm) WRITE ( numond, namtra_adv )
179
180      IF(lwp) THEN                    ! Namelist print
181         WRITE(numout,*)
182         WRITE(numout,*) 'tra_adv_init : choice/control of the tracer advection scheme'
183         WRITE(numout,*) '~~~~~~~~~~~'
184         WRITE(numout,*) '   Namelist namtra_adv : chose a advection scheme for tracers'
185         WRITE(numout,*) '      2nd order advection scheme     ln_traadv_cen2    = ', ln_traadv_cen2
186         WRITE(numout,*) '      TVD advection scheme           ln_traadv_tvd     = ', ln_traadv_tvd
187         WRITE(numout,*) '      MUSCL  advection scheme        ln_traadv_muscl   = ', ln_traadv_muscl
188         WRITE(numout,*) '      MUSCL2 advection scheme        ln_traadv_muscl2  = ', ln_traadv_muscl2
189         WRITE(numout,*) '      UBS    advection scheme        ln_traadv_ubs     = ', ln_traadv_ubs
190         WRITE(numout,*) '      QUICKEST advection scheme      ln_traadv_qck     = ', ln_traadv_qck
191         WRITE(numout,*) '      upstream scheme within muscl   ln_traadv_msc_ups = ', ln_traadv_msc_ups
192      ENDIF
193
194      ioptio = 0                      ! Parameter control
195      IF( ln_traadv_cen2   )   ioptio = ioptio + 1
196      IF( ln_traadv_tvd    )   ioptio = ioptio + 1
197      IF( ln_traadv_muscl  )   ioptio = ioptio + 1
198      IF( ln_traadv_muscl2 )   ioptio = ioptio + 1
199      IF( ln_traadv_ubs    )   ioptio = ioptio + 1
200      IF( ln_traadv_qck    )   ioptio = ioptio + 1
201      IF( lk_esopa         )   ioptio =          1
202
203      IF( ioptio /= 1 )   CALL ctl_stop( 'Choose ONE advection scheme in namelist namtra_adv' )
204
205      !                              ! Set nadv
206      IF( ln_traadv_cen2   )   nadv =  1
207      IF( ln_traadv_tvd    )   nadv =  2
208      IF( ln_traadv_muscl  )   nadv =  3
209      IF( ln_traadv_muscl2 )   nadv =  4
210      IF( ln_traadv_ubs    )   nadv =  5
211      IF( ln_traadv_qck    )   nadv =  6
212      IF( lk_esopa         )   nadv = -1
213
214      IF(lwp) THEN                   ! Print the choice
215         WRITE(numout,*)
216         IF( nadv ==  1 )   WRITE(numout,*) '         2nd order scheme is used'
217         IF( nadv ==  2 )   WRITE(numout,*) '         TVD       scheme is used'
218         IF( nadv ==  3 )   WRITE(numout,*) '         MUSCL     scheme is used'
219         IF( nadv ==  4 )   WRITE(numout,*) '         MUSCL2    scheme is used'
220         IF( nadv ==  5 )   WRITE(numout,*) '         UBS       scheme is used'
221         IF( nadv ==  6 )   WRITE(numout,*) '         QUICKEST  scheme is used'
222         IF( nadv == -1 )   WRITE(numout,*) '         esopa test: use all advection scheme'
223      ENDIF
224      !
225      CALL tra_adv_mle_init          ! initialisation of the Mixed Layer Eddy parametrisation (MLE)
226      !
227   END SUBROUTINE tra_adv_init
228
229  !!======================================================================
230END MODULE traadv
Note: See TracBrowser for help on using the repository browser.