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.
dyndmp.F90 in NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/C1D – NEMO

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/C1D/dyndmp.F90 @ 12340

Last change on this file since 12340 was 12340, checked in by acc, 4 years ago

Branch 2019/dev_r11943_MERGE_2019. This commit introduces basic do loop macro
substitution to the 2019 option 1, merge branch. These changes have been SETTE
tested. The only addition is the do_loop_substitute.h90 file in the OCE directory but
the macros defined therein are used throughout the code to replace identifiable, 2D-
and 3D- nested loop opening and closing statements with single-line alternatives. Code
indents are also adjusted accordingly.

The following explanation is taken from comments in the new header file:

This header file contains preprocessor definitions and macros used in the do-loop
substitutions introduced between version 4.0 and 4.2. The primary aim of these macros
is to assist in future applications of tiling to improve performance. This is expected
to be achieved by alternative versions of these macros in selected locations. The
initial introduction of these macros simply replaces all identifiable nested 2D- and
3D-loops with single line statements (and adjusts indenting accordingly). Do loops
are identifiable if they comform to either:

DO jk = ....

DO jj = .... DO jj = ...

DO ji = .... DO ji = ...
. OR .
. .

END DO END DO

END DO END DO

END DO

and white-space variants thereof.

Additionally, only loops with recognised jj and ji loops limits are treated; these are:
Lower limits of 1, 2 or fs_2
Upper limits of jpi, jpim1 or fs_jpim1 (for ji) or jpj, jpjm1 or fs_jpjm1 (for jj)

The macro naming convention takes the form: DO_2D_BT_LR where:

B is the Bottom offset from the PE's inner domain;
T is the Top offset from the PE's inner domain;
L is the Left offset from the PE's inner domain;
R is the Right offset from the PE's inner domain

So, given an inner domain of 2,jpim1 and 2,jpjm1, a typical example would replace:

DO jj = 2, jpj

DO ji = 1, jpim1
.
.

END DO

END DO

with:

DO_2D_01_10
.
.
END_2D

similar conventions apply to the 3D loops macros. jk loop limits are retained
through macro arguments and are not restricted. This includes the possibility of
strides for which an extra set of DO_3DS macros are defined.

In the example definition below the inner PE domain is defined by start indices of
(kIs, kJs) and end indices of (kIe, KJe)

#define DO_2D_00_00 DO jj = kJs, kJe ; DO ji = kIs, kIe
#define END_2D END DO ; END DO

TO DO:


Only conventional nested loops have been identified and replaced by this step. There are constructs such as:

DO jk = 2, jpkm1

z2d(:,:) = z2d(:,:) + e3w(:,:,jk,Kmm) * z3d(:,:,jk) * wmask(:,:,jk)

END DO

which may need to be considered.

  • Property svn:keywords set to Id
File size: 10.6 KB
Line 
1MODULE dyndmp
2   !!======================================================================
3   !!                       ***  MODULE  dyndmp  ***
4   !! Ocean dynamics: internal restoring trend on momentum (U and V current)
5   !!                 This should only be used for C1D case in current form
6   !!======================================================================
7   !! History :  3.5  ! 2013-08  (D. Calvert)  Original code
8   !!            3.6  ! 2014-08  (T. Graham) Modified to use netcdf file of
9   !!                             restoration coefficients supplied to tradmp
10   !!----------------------------------------------------------------------
11
12   !!----------------------------------------------------------------------
13   !!   dyn_dmp_alloc : allocate dyndmp arrays
14   !!   dyn_dmp_init  : namelist read, parameter control and resto coeff.
15   !!   dyn_dmp       : update the momentum trend with the internal damping
16   !!----------------------------------------------------------------------
17   USE oce            ! ocean: variables
18   USE dom_oce        ! ocean: domain variables
19   USE c1d            ! 1D vertical configuration
20   USE tradmp         ! ocean: internal damping
21   USE zdf_oce        ! ocean: vertical physics
22   USE phycst         ! physical constants
23   USE dtauvd         ! data: U & V current
24   USE zdfmxl         ! vertical physics: mixed layer depth
25   !
26   USE in_out_manager ! I/O manager
27   USE lib_mpp        ! MPP library
28   USE prtctl         ! Print control
29   USE timing         ! Timing
30   USE iom            ! I/O manager
31
32   IMPLICIT NONE
33   PRIVATE
34
35   PUBLIC   dyn_dmp_init ! routine called by nemogcm.F90
36   PUBLIC   dyn_dmp      ! routine called by step_c1d.F90
37
38   LOGICAL, PUBLIC ::   ln_dyndmp   !: Flag for Newtonian damping
39
40   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::  utrdmp    !: damping U current trend (m/s2)
41   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::  vtrdmp    !: damping V current trend (m/s2)
42   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::  resto_uv  !: restoring coeff. on U & V current
43
44   !! * Substitutions
45#  include "vectopt_loop_substitute.h90"
46#  include "do_loop_substitute.h90"
47   !!----------------------------------------------------------------------
48   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
49   !! $Id$
50   !! Software governed by the CeCILL license (see ./LICENSE)
51   !!----------------------------------------------------------------------
52CONTAINS
53
54   INTEGER FUNCTION dyn_dmp_alloc()
55      !!----------------------------------------------------------------------
56      !!                ***  FUNCTION dyn_dmp_alloc  ***
57      !!----------------------------------------------------------------------
58      ALLOCATE( utrdmp(jpi,jpj,jpk), vtrdmp(jpi,jpj,jpk), resto_uv(jpi,jpj,jpk), STAT= dyn_dmp_alloc )
59      !
60      CALL mpp_sum ( 'dyndmp', dyn_dmp_alloc )
61      IF( dyn_dmp_alloc > 0 )   CALL ctl_warn('dyn_dmp_alloc: allocation of arrays failed')
62      !
63   END FUNCTION dyn_dmp_alloc
64
65
66   SUBROUTINE dyn_dmp_init
67      !!----------------------------------------------------------------------
68      !!                  ***  ROUTINE dyn_dmp_init  ***
69      !!
70      !! ** Purpose :   Initialization for the Newtonian damping
71      !!
72      !! ** Method  : - read the ln_dyndmp parameter from the namc1d_dyndmp namelist
73      !!              - allocate damping arrays
74      !!              - check the parameters of the namtra_dmp namelist
75      !!              - calculate damping coefficient
76      !!----------------------------------------------------------------------
77      INTEGER ::   ios, imask   ! local integers
78      !!
79      NAMELIST/namc1d_dyndmp/ ln_dyndmp
80      !!----------------------------------------------------------------------
81      !
82      READ  ( numnam_ref, namc1d_dyndmp, IOSTAT = ios, ERR = 901)
83901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namc1d_dyndmp in reference namelist' )
84      READ  ( numnam_cfg, namc1d_dyndmp, IOSTAT = ios, ERR = 902 )
85902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namc1d_dyndmp in configuration namelist' )
86      IF(lwm) WRITE ( numond, namc1d_dyndmp )
87      !
88      IF(lwp) THEN                           ! control print
89         WRITE(numout,*)
90         WRITE(numout,*) 'dyn_dmp_init : U and V current Newtonian damping'
91         WRITE(numout,*) '~~~~~~~~~~~~'
92         WRITE(numout,*) '   Namelist namc1d_dyndmp : Set damping flag'
93         WRITE(numout,*) '      add a damping term or not       ln_dyndmp = ', ln_dyndmp
94         WRITE(numout,*) '   Namelist namtra_dmp    : Set damping parameters'
95         WRITE(numout,*) '      Apply relaxation   or not       ln_tradmp = ', ln_tradmp
96         WRITE(numout,*) '      mixed layer damping option      nn_zdmp   = ', nn_zdmp
97         WRITE(numout,*) '      Damping file name               cn_resto  = ', cn_resto
98         WRITE(numout,*)
99      ENDIF
100      !
101      IF( ln_dyndmp ) THEN
102         !                                   !==   allocate the data arrays   ==!
103         IF( dyn_dmp_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dyn_dmp_init: unable to allocate arrays' )
104         !
105         SELECT CASE ( nn_zdmp )             !==   control print of vertical option   ==!
106         CASE ( 0    )   ;   IF(lwp) WRITE(numout,*) '   momentum damping throughout the water column'
107         CASE ( 1    )   ;   IF(lwp) WRITE(numout,*) '   no momentum damping in the turbocline (avt > 5 cm2/s)'
108         CASE ( 2    )   ;   IF(lwp) WRITE(numout,*) '   no momentum damping in the mixed layer'
109         CASE DEFAULT
110            WRITE(ctmp1,*) '          bad flag value for nn_zdmp = ', nn_zdmp
111            CALL ctl_stop(ctmp1)
112         END SELECT
113         !
114         IF( .NOT. ln_uvd_dyndmp ) THEN      ! force the initialization of U & V current data for damping
115            CALL ctl_warn( 'dyn_dmp_init: U & V current read data not initialized, we force ln_uvd_dyndmp=T' )
116            CALL dta_uvd_init( ld_dyndmp=ln_dyndmp )
117         ENDIF
118         !
119         utrdmp(:,:,:) = 0._wp               ! internal damping trends
120         vtrdmp(:,:,:) = 0._wp
121         !
122         !Read in mask from file
123         CALL iom_open ( cn_resto, imask)
124         CALL iom_get  ( imask, jpdom_autoglo, 'resto', resto)
125         CALL iom_close( imask )
126      ENDIF
127      !
128   END SUBROUTINE dyn_dmp_init
129
130
131   SUBROUTINE dyn_dmp( kt, Kbb, Kmm, puu, pvv, Krhs )
132      !!----------------------------------------------------------------------
133      !!                   ***  ROUTINE dyn_dmp  ***
134      !!                 
135      !! ** Purpose :   Compute the momentum trends due to a newtonian damping
136      !!      of the ocean velocities towards the given data and add it to the
137      !!      general momentum trends.
138      !!
139      !! ** Method  :   Compute Newtonian damping towards u_dta and v_dta
140      !!      and add to the general momentum trends:
141      !!                     puu(Krhs) = puu(Krhs) + resto_uv * (u_dta - puu(Kbb))
142      !!                     pvv(Krhs) = pvv(Krhs) + resto_uv * (v_dta - pvv(Kbb))
143      !!      The trend is computed either throughout the water column
144      !!      (nn_zdmp=0), where the vertical mixing is weak (nn_zdmp=1) or
145      !!      below the well mixed layer (nn_zdmp=2)
146      !!
147      !! ** Action  : - (puu(:,:,:,Krhs),pvv(:,:,:,Krhs))   momentum trends updated with the damping trend
148      !!----------------------------------------------------------------------
149      INTEGER                             , INTENT(in   ) ::   kt             ! ocean time-step index
150      INTEGER                             , INTENT(in   ) ::   Kbb, Kmm, Krhs ! ocean time level indices
151      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) ::   puu, pvv       ! ocean velocities and RHS of momentum equation
152      !!
153      INTEGER  ::   ji, jj, jk   ! dummy loop indices
154      REAL(wp) ::   zua, zva     ! local scalars
155      REAL(wp), DIMENSION(jpi,jpj,jpk,2) ::   zuv_dta   ! Read in data
156      !!----------------------------------------------------------------------
157      !
158      IF( ln_timing )   CALL timing_start( 'dyn_dmp' )
159      !
160      !
161      !                           !==   read and interpolate U & V current data at kt   ==!
162      CALL dta_uvd( kt, Kmm, zuv_dta ) !!! NOTE: This subroutine must be altered for use outside
163                                  !!!       the C1D context (use of U,V grid variables)
164      !
165      SELECT CASE ( nn_zdmp )     !==   Calculate/add Newtonian damping to the momentum trend   ==!
166      !
167      CASE( 0 )                   ! Newtonian damping throughout the water column
168         DO_3D_00_00( 1, jpkm1 )
169            zua = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,1) - puu(ji,jj,jk,Kbb) )
170            zva = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,2) - pvv(ji,jj,jk,Kbb) )
171            puu(ji,jj,jk,Krhs) = puu(ji,jj,jk,Krhs) + zua
172            pvv(ji,jj,jk,Krhs) = pvv(ji,jj,jk,Krhs) + zva
173            utrdmp(ji,jj,jk) = zua           ! save the trends
174            vtrdmp(ji,jj,jk) = zva     
175         END_3D
176         !
177      CASE ( 1 )                  ! no damping above the turbocline (avt > 5 cm2/s)
178         DO_3D_00_00( 1, jpkm1 )
179            IF( avt(ji,jj,jk) <= avt_c ) THEN
180               zua = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,1) - puu(ji,jj,jk,Kbb) )
181               zva = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,2) - pvv(ji,jj,jk,Kbb) )
182            ELSE
183               zua = 0._wp
184               zva = 0._wp 
185            ENDIF
186            puu(ji,jj,jk,Krhs) = puu(ji,jj,jk,Krhs) + zua
187            pvv(ji,jj,jk,Krhs) = pvv(ji,jj,jk,Krhs) + zva
188            utrdmp(ji,jj,jk) = zua           ! save the trends
189            vtrdmp(ji,jj,jk) = zva
190         END_3D
191         !
192      CASE ( 2 )                  ! no damping in the mixed layer
193         DO_3D_00_00( 1, jpkm1 )
194            IF( gdept(ji,jj,jk,Kmm) >= hmlp (ji,jj) ) THEN
195               zua = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,1) - puu(ji,jj,jk,Kbb) )
196               zva = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,2) - pvv(ji,jj,jk,Kbb) )
197            ELSE
198               zua = 0._wp
199               zva = 0._wp 
200            ENDIF
201            puu(ji,jj,jk,Krhs) = puu(ji,jj,jk,Krhs) + zua
202            pvv(ji,jj,jk,Krhs) = pvv(ji,jj,jk,Krhs) + zva
203            utrdmp(ji,jj,jk) = zua           ! save the trends
204            vtrdmp(ji,jj,jk) = zva
205         END_3D
206         !
207      END SELECT
208      !
209      !                           ! Control print
210      IF( sn_cfctl%l_prtctl   )   CALL prt_ctl( tab3d_1=puu(:,:,:,Krhs), clinfo1=' dmp  - Ua: ', mask1=umask,   &
211         &                                      tab3d_2=pvv(:,:,:,Krhs), clinfo2=       ' Va: ', mask2=vmask, clinfo3='dyn' )
212      !
213      !
214      IF( ln_timing )   CALL timing_stop( 'dyn_dmp')
215      !
216   END SUBROUTINE dyn_dmp
217
218   !!======================================================================
219END MODULE dyndmp
Note: See TracBrowser for help on using the repository browser.