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 @ 14053

Last change on this file since 14053 was 14053, checked in by techene, 4 years ago

#2385 added to the trunk

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