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.
usrdef_istate.F90 in NEMO/trunk/tests/CANAL/MY_SRC – NEMO

source: NEMO/trunk/tests/CANAL/MY_SRC/usrdef_istate.F90 @ 14433

Last change on this file since 14433 was 14433, checked in by smasson, 3 years ago

trunk: merge dev_r14312_MPI_Interface into the trunk, #2598

  • Property svn:keywords set to Id
File size: 15.6 KB
Line 
1MODULE usrdef_istate
2   !!======================================================================
3   !!                     ***  MODULE usrdef_istate   ***
4   !!
5   !!                      ===  CANAL configuration  ===
6   !!
7   !! User defined : set the initial state of a user configuration
8   !!======================================================================
9   !! History :  NEMO ! 2017-11  (J. Chanut) Original code
10   !!----------------------------------------------------------------------
11
12   !!----------------------------------------------------------------------
13   !!  usr_def_istate : initial state in Temperature and salinity
14   !!----------------------------------------------------------------------
15   USE par_oce        ! ocean space and time domain
16   USE dom_oce       
17   USE phycst         ! physical constants
18   USE eosbn2, ONLY: rn_a0
19   !
20   USE in_out_manager ! I/O manager
21   USE lib_mpp        ! MPP library
22   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
23   !   
24   USE usrdef_nam
25   
26   IMPLICIT NONE
27   PRIVATE
28
29   PUBLIC   usr_def_istate       ! called by istate.F90
30   PUBLIC   usr_def_istate_ssh   ! called by sshwzv.F90
31
32   !! * Substitutions
33#  include "do_loop_substitute.h90"
34   !!----------------------------------------------------------------------
35   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
36   !! $Id$
37   !! Software governed by the CeCILL license (see ./LICENSE)
38   !!----------------------------------------------------------------------
39CONTAINS
40 
41   SUBROUTINE usr_def_istate( pdept, ptmask, pts, pu, pv )
42      !!----------------------------------------------------------------------
43      !!                   ***  ROUTINE usr_def_istate  ***
44      !!
45      !! ** Purpose :   Initialization of the dynamics and tracers
46      !!                Here CANAL configuration
47      !!
48      !! ** Method  :   Set a gaussian anomaly of pressure and associated
49      !!                geostrophic velocities
50      !!----------------------------------------------------------------------
51      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(in   ) ::   pdept   ! depth of t-point               [m]
52      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(in   ) ::   ptmask  ! t-point ocean mask             [m]
53      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(  out) ::   pts     ! T & S fields      [Celsius ; g/kg]
54      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(  out) ::   pu      ! i-component of the velocity  [m/s]
55      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(  out) ::   pv      ! j-component of the velocity  [m/s]
56      !
57      INTEGER  :: ji, jj, jk, jl  ! dummy loop indices
58      REAL(wp) :: zx, zy, zP0, zumax, zlambda, zr_lambda2, zn2, zf0, zH, zrho1, za, zf, zdzF
59      REAL(wp) :: zpsurf, zdyPs, zdxPs
60      REAL(wp) :: zdt, zdu, zdv
61      REAL(wp) :: zjetx, zjety, zbeta
62      REAL(wp), DIMENSION(jpi,jpj)  ::   zrandom
63      !!----------------------------------------------------------------------
64      !
65      IF(lwp) WRITE(numout,*)
66      IF(lwp) WRITE(numout,*) 'usr_def_istate : CANAL configuration, analytical definition of initial state'
67      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~   '
68      !
69      zjetx = ABS(rn_ujetszx)/2.
70      zjety = ABS(rn_ujetszy)/2.
71      !
72      zf0   = 2._wp * omega * SIN( rad * rn_ppgphi0 )
73      !
74      SELECT CASE(nn_initcase)
75
76      CASE(-1)    ! stratif at rest
77
78         ! temperature:
79         pts(:,:,1,jp_tem) = 25. !!30._wp
80         pts(:,:,2:jpk,jp_tem) = 22. !!24._wp
81         ! salinity: 
82         pts(:,:,:,jp_sal) = 35._wp
83         ! velocities:
84         pu(:,:,:) = 0.
85         pv(:,:,:) = 0.
86
87      CASE(0)    ! rest
88         !
89         ! temperature:
90         pts(:,:,:,jp_tem) = 10._wp
91         ! salinity: 
92         pts(:,:,:,jp_sal) = 35._wp
93         ! velocities:
94         pu(:,:,:) = 0.
95         pv(:,:,:) = 0.
96         
97      CASE(1)    ! geostrophic zonal jet from -zjety to +zjety
98         !
99         ! temperature:
100         pts(:,:,:,jp_tem) = 10._wp
101         ! salinity: 
102         pts(:,:,jpk,jp_sal) = 0.
103         DO jk=1, jpkm1
104            WHERE( ABS(gphit) <= zjety )
105!!$            WHERE( ABS(gphit) <= zjety*0.5 .AND. ABS(glamt) <= zjety*0.5 ) ! for a square of salt
106               pts(:,:,jk,jp_sal) = 35.
107            ELSEWHERE
108               pts(:,:,jk,jp_sal) = 30.
109            END WHERE                   
110         END DO
111         ! velocities:
112         pu(:,:,:) = 0.
113         DO jk=1, jpkm1
114            WHERE( ABS(gphit) <= zjety ) pu(:,:,jk) = rn_uzonal
115         END DO
116         pv(:,:,:) = 0.
117         !                 
118      CASE(2)    ! geostrophic zonal current shear
119         !
120         ! temperature:
121         pts(:,:,:,jp_tem) = 10._wp
122         ! salinity: 
123         pts(:,:,:,jp_sal) = 30.
124         DO jk=1, jpkm1
125            WHERE( ABS(gphiv) <= zjety ) pts(:,:,jk,jp_sal) = 30. + SIGN(1.,gphiv(:,:))
126         END DO
127         ! velocities:
128         pu(:,:,:) = 0.
129         DO jk=1, jpkm1
130            WHERE( ABS(gphiv) <= zjety ) pu(:,:,jk) = SIGN(rn_uzonal,gphit(:,:))*SIGN(1.,rn_uzonal)
131            WHERE( ABS(gphiv) == 0.    ) pu(:,:,jk) = 0. 
132         END DO
133         pv(:,:,:) = 0.
134         !                 
135      CASE(3)    ! gaussian zonal currant
136         !
137         ! zonal current
138         DO jk=1, jpkm1
139            ! gphit and lambda are both in km
140            pu(:,:,jk) = rn_uzonal * EXP( - 0.5 * gphit(:,:)**2 / rn_lambda**2 )
141         END DO
142         ! temperature:
143         pts(:,:,:,jp_tem) = 10._wp
144         ! salinity: 
145         DO jk=1, jpkm1
146            pts(:,:,jk,jp_sal) = gphit(:,:)
147         END DO
148         ! velocities:
149         pv(:,:,:) = 0.
150         !           
151      CASE(4)    ! geostrophic zonal pulse
152         !
153         DO_2D( 1, 1, 1, 1 )
154            IF ( ABS(glamt(ji,jj)) <= zjetx ) THEN
155               zdu = rn_uzonal
156            ELSEIF ( ABS(glamt(ji,jj)) <= zjetx + 100. ) THEN
157               zdu = rn_uzonal * ( ( zjetx-ABS(glamt(ji,jj)) )/100. + 1. )
158            ELSE
159               zdu = 0.
160            ENDIF
161            IF ( ABS(gphit(ji,jj)) <= zjety ) THEN
162               pu(ji,jj,:) = zdu
163               pts(ji,jj,:,jp_sal) = zdu / rn_uzonal + 1.
164            ELSE
165               pu(ji,jj,:) = 0.
166               pts(ji,jj,:,jp_sal) = 1.
167            ENDIF
168         END_2D
169         !
170         ! temperature:
171         pts(:,:,:,jp_tem) = 10._wp * ptmask(:,:,:)       
172         pv(:,:,:) = 0.
173         !
174      CASE(5)    ! vortex
175         !
176         zf0   = 2._wp * omega * SIN( rad * rn_ppgphi0 )
177         zumax = rn_vtxmax * SIGN(1._wp, zf0)  ! Here Anticyclonic: set zumax=-1 for cyclonic
178         zlambda = SQRT(2._wp)*rn_lambda*1.e3       ! Horizontal scale in meters
179         zn2 = 3.e-3**2
180         zH = 0.5_wp * 5000._wp
181         !
182         zr_lambda2 = 1._wp / zlambda**2
183         zP0 = rho0 * zf0 * zumax * zlambda * SQRT(EXP(1._wp)/2._wp)
184         !
185         DO_2D( 1, 1, 1, 1 )
186            zx = glamt(ji,jj) * 1.e3
187            zy = gphit(ji,jj) * 1.e3
188            ! Surface pressure: P(x,y,z) = F(z) * Psurf(x,y)
189            zpsurf = zP0 * EXP(-(zx**2+zy**2)*zr_lambda2) - rho0 * ff_t(ji,jj) * rn_uzonal * zy
190            ! temperature:
191            DO jk=1,jpk
192               zdt =  pdept(ji,jj,jk) 
193               zrho1 = rho0 * (1._wp + zn2*zdt/grav)
194               IF (zdt < zH) THEN
195                  zdzF = (1._wp-EXP(zdt-zH)) / (zH-1._wp + EXP(-zH))   ! F'(z)
196                  zrho1 = zrho1 - zdzF * zpsurf / grav    ! -1/g Dz(P) = -1/g * F'(z) * Psurf(x,y)
197               ENDIF
198               !               pts(ji,jj,jk,jp_tem) = (20._wp + (rho0-zrho1) / rn_a0) * ptmask(ji,jj,jk)
199               pts(ji,jj,jk,jp_tem) = (10._wp + (rho0-zrho1) / rn_a0) * ptmask(ji,jj,jk)
200            END DO
201         END_2D
202         !
203         ! salinity: 
204         pts(:,:,:,jp_sal) = 35._wp * ptmask(:,:,:) 
205         !
206         ! velocities:
207         za = 2._wp * zP0 / zlambda**2
208         DO_2D( 0, 0, 0, 0 )
209            zx = glamu(ji,jj) * 1.e3
210            zy = gphiu(ji,jj) * 1.e3
211            DO jk=1, jpk
212               zdu = 0.5_wp * (pdept(ji,jj,jk) + pdept(ji+1,jj,jk))
213               IF (zdu < zH) THEN
214                  zf = (zH-1._wp-zdu+EXP(zdu-zH)) / (zH-1._wp+EXP(-zH))
215                  zdyPs = - za * zy * EXP(-(zx**2+zy**2)*zr_lambda2) - rho0 * ff_t(ji,jj) * rn_uzonal
216                  pu(ji,jj,jk) = - zf / ( rho0 * ff_t(ji,jj) ) * zdyPs * ptmask(ji,jj,jk) * ptmask(ji+1,jj,jk)
217               ELSE
218                  pu(ji,jj,jk) = 0._wp
219               ENDIF
220            END DO
221         END_2D
222         !
223         DO_2D( 0, 0, 0, 0 )
224            zx = glamv(ji,jj) * 1.e3
225            zy = gphiv(ji,jj) * 1.e3
226            DO jk=1, jpk
227               zdv = 0.5_wp * (pdept(ji,jj,jk) + pdept(ji,jj+1,jk))
228               IF (zdv < zH) THEN
229                  zf = (zH-1._wp-zdv+EXP(zdv-zH)) / (zH-1._wp+EXP(-zH))
230                  zdxPs = - za * zx * EXP(-(zx**2+zy**2)*zr_lambda2)
231                  pv(ji,jj,jk) = zf / ( rho0 * ff_f(ji,jj) ) * zdxPs * ptmask(ji,jj,jk) * ptmask(ji,jj+1,jk)
232               ELSE
233                  pv(ji,jj,jk) = 0._wp
234               ENDIF
235            END DO
236         END_2D
237         !           
238      END SELECT
239      !
240      CALL lbc_lnk( 'usrdef_istate', pts , 'T',  1. )
241      CALL lbc_lnk( 'usrdef_istate', pu, 'U', -1., pv, 'V', -1. )
242
243   END SUBROUTINE usr_def_istate
244
245 
246   SUBROUTINE usr_def_istate_ssh( ptmask, pssh )
247      !!----------------------------------------------------------------------
248      !!                   ***  ROUTINE usr_def_istate_ssh  ***
249      !!
250      !! ** Purpose :   Initialization of the dynamics and tracers
251      !!                Here CANAL configuration
252      !!
253      !! ** Method  :   Set ssh
254      !!----------------------------------------------------------------------
255      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(in   ) ::   ptmask  ! t-point ocean mask             [m]
256      REAL(wp), DIMENSION(jpi,jpj)         , INTENT(  out) ::   pssh    ! sea-surface height
257      !
258      INTEGER  :: ji, jj, jk, jl  ! dummy loop indices
259      REAL(wp) :: zx, zy, zP0, zumax, zlambda, zr_lambda2, zn2, zf0, zH, zrho1, za, zf, zdzF
260      REAL(wp) :: zpsurf, zdyPs, zdxPs
261      REAL(wp) :: zdt, zdu, zdv
262      REAL(wp) :: zjetx, zjety, zbeta
263      REAL(wp), DIMENSION(jpi,jpj)  ::   zrandom
264      !!----------------------------------------------------------------------
265      !
266      IF(lwp) WRITE(numout,*)
267      IF(lwp) WRITE(numout,*) 'usr_def_istate_ssh : CANAL configuration, analytical definition of initial state'
268      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~   '
269      !
270      IF (ln_sshnoise) CALL RANDOM_NUMBER(zrandom)
271      zjetx = ABS(rn_ujetszx)/2.
272      zjety = ABS(rn_ujetszy)/2.
273      !
274      SELECT CASE(nn_initcase)
275      CASE(0)                      !==   rest  ==!
276         !
277         pssh(:,:) = 0.
278         !
279      CASE(1)                      !==  geostrophic zonal jet from -zjety to +zjety  ==!
280         !
281         SELECT CASE( nn_fcase )
282         CASE(0)                          !* f = f0 : ssh = - fuy / g
283            WHERE( ABS(gphit) <= zjety )
284               pssh(:,:) = - ff_t(:,:) * rn_uzonal * gphit(:,:) * 1.e3 / grav
285            ELSEWHERE
286               pssh(:,:) = - ff_t(:,:) * rn_uzonal * SIGN(zjety, gphit(:,:)) * 1.e3 / grav
287            END WHERE
288         CASE(1)                          !* f = f0 + beta*y : ssh = - u / g * ( fy + 0.5 * beta * y^2 )
289            zbeta = 2._wp * omega * COS( rad * rn_ppgphi0 ) / ra
290            WHERE( ABS(gphit) <= zjety )
291               pssh(:,:) = - rn_uzonal / grav * ( ff_t(:,:) * gphit(:,:) * 1.e3 + 0.5 * zbeta * gphit(:,:) * gphit(:,:) * 1.e6 )
292            ELSEWHERE
293               pssh(:,:) = - rn_uzonal / grav * ( ff_t(:,:) * SIGN(zjety, gphit(:,:)) * 1.e3   &
294                  &                             + 0.5 * zbeta * zjety * zjety * 1.e6 )
295            END WHERE
296         END SELECT
297         !                 
298      CASE(2)                      !==  geostrophic zonal current shear  ==!
299         !
300         SELECT CASE( nn_fcase )
301         CASE(0)                          !* f = f0 : ssh = - fuy / g
302            WHERE( ABS(gphit) <= zjety )
303               pssh(:,:) = - ff_t(:,:) * rn_uzonal * ABS(gphit(:,:)) * 1.e3 / grav
304            ELSEWHERE
305               pssh(:,:) = - ff_t(:,:) * rn_uzonal * zjety * 1.e3 / grav
306            END WHERE
307         CASE(1)                          !* f = f0 + beta*y : ssh = - u / g * ( fy + 0.5 * beta * y^2 )
308            zbeta = 2._wp * omega * COS( rad * rn_ppgphi0 ) / ra
309            WHERE( ABS(gphit) <= zjety )
310               pssh(:,:) = - SIGN(rn_uzonal, gphit(:,:)) / grav   &
311                  &        * ( ff_t(:,:) * gphit(:,:) * 1.e3 + 0.5 * zbeta * gphit(:,:) * gphit(:,:) * 1.e6 )
312            ELSEWHERE
313               pssh(:,:) = - SIGN(rn_uzonal, gphit(:,:)) / grav   &
314                  &        * ( ff_t(:,:) * SIGN(zjety, gphit(:,:)) * 1.e3 + 0.5 * zbeta * zjety * zjety * 1.e6 )
315            END WHERE
316         END SELECT
317         !                 
318      CASE(3)                      !==  gaussian zonal currant  ==!
319         !
320         pssh(:,1) = - ff_t(:,1) / grav * e2t(:,1) * rn_uzonal * EXP( - 0.5 * gphit(:,1)**2 / rn_lambda**2 )
321         DO jl=1, jpnj
322            DO_2D( 0, 0, 0, 0 )
323               pssh(ji,jj) = pssh(ji,jj-1) - ff_t(ji,jj) / grav * rn_uzonal * EXP( - 0.5 * gphit(ji,jj)**2 / rn_lambda**2 ) * e2t(ji,jj)
324            END_2D
325            CALL lbc_lnk( 'usrdef_istate_ssh', pssh, 'T',  1. )
326         END DO
327         !           
328      CASE(4)                      !==  geostrophic zonal pulse !!st need to implement a way to separate ssh properly  ==!
329         !
330         DO_2D( 1, 1, 1, 1 )
331            IF ( ABS(glamt(ji,jj)) <= zjetx ) THEN
332               zdu = rn_uzonal
333            ELSEIF ( ABS(glamt(ji,jj)) <= zjetx + 100. ) THEN
334               zdu = rn_uzonal * ( ( zjetx-ABS(glamt(ji,jj)) )/100. + 1. )
335            ELSE
336               zdu = 0.
337            ENDIF
338            IF ( ABS(gphit(ji,jj)) <= zjety ) THEN
339               pssh(ji,jj) = - ff_t(ji,jj) * zdu * gphit(ji,jj) * 1.e3 / grav
340            ELSE
341               pssh(ji,jj) = - ff_t(ji,jj) * zdu * SIGN(zjety,gphit(ji,jj)) * 1.e3 / grav 
342            ENDIF
343         END_2D
344         !
345      CASE(5)                    !==  vortex  ==!
346         !
347         zf0   = 2._wp * omega * SIN( rad * rn_ppgphi0 )
348         zumax = rn_vtxmax * SIGN(1._wp, zf0)   ! Here Anticyclonic: set zumax=-1 for cyclonic
349         zlambda = SQRT(2._wp)*rn_lambda        ! Horizontal scale in meters
350         zn2 = 3.e-3**2
351         zH = 0.5_wp * 5000._wp
352         !
353         zr_lambda2 = 1._wp / zlambda**2
354         zP0 = rho0 * zf0 * zumax * zlambda * SQRT(EXP(1._wp)/2._wp)
355         !
356         DO_2D( 1, 1, 1, 1 )
357            zx = glamt(ji,jj) * 1.e3
358            zy = gphit(ji,jj) * 1.e3
359            !                                   ! Surface pressure: P(x,y,z) = F(z) * Psurf(x,y)
360            zpsurf = zP0 * EXP(-(zx**2+zy**2)*zr_lambda2) - rho0 * ff_t(ji,jj) * rn_uzonal * zy
361            pssh(ji,jj) = 0.
362            DO jl=1,5
363               zdt = pssh(ji,jj)
364               zdzF = (1._wp - EXP(zdt-zH)) / (zH - 1._wp + EXP(-zH))   ! F'(z)
365               zrho1 = rho0 * (1._wp + zn2*zdt/grav) - zdzF * zpsurf / grav    ! -1/g Dz(P) = -1/g * F'(z) * Psurf(x,y)
366               pssh(ji,jj) = zpsurf / (zrho1*grav) * ptmask(ji,jj,1)   ! ssh = Psurf / (Rho*g)
367            END DO
368         END_2D
369         !           
370      END SELECT
371      !                          !==  add noise  ==!
372      IF (ln_sshnoise) THEN
373         CALL RANDOM_SEED()
374         CALL RANDOM_NUMBER(zrandom)
375         pssh(:,:) = pssh(:,:) + ( 0.1  * zrandom(:,:) - 0.05 )
376      ENDIF
377      CALL lbc_lnk( 'usrdef_istate_ssh', pssh, 'T',  1. )
378      !
379   END SUBROUTINE usr_def_istate_ssh
380   
381   !!======================================================================
382END MODULE usrdef_istate
Note: See TracBrowser for help on using the repository browser.