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/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/tests/CANAL/MY_SRC – NEMO

source: NEMO/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/tests/CANAL/MY_SRC/usrdef_istate.F90 @ 14200

Last change on this file since 14200 was 14200, checked in by mcastril, 3 years ago

Merging r14117 through r14199 into dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final

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