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/r12377_ticket2386/tests/CANAL/MY_SRC – NEMO

source: NEMO/branches/2020/r12377_ticket2386/tests/CANAL/MY_SRC/usrdef_istate.F90 @ 13540

Last change on this file since 13540 was 13540, checked in by andmirek, 4 years ago

Ticket #2386: update to latest trunk

  • Property svn:keywords set to Id
File size: 12.4 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
30   !! * Substitutions
31#  include "do_loop_substitute.h90"
32   !!----------------------------------------------------------------------
33   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
34   !! $Id$
35   !! Software governed by the CeCILL license (see ./LICENSE)
36   !!----------------------------------------------------------------------
37CONTAINS
38 
39   SUBROUTINE usr_def_istate( pdept, ptmask, pts, pu, pv, pssh )
40      !!----------------------------------------------------------------------
41      !!                   ***  ROUTINE usr_def_istate  ***
42      !!
43      !! ** Purpose :   Initialization of the dynamics and tracers
44      !!                Here CANAL configuration
45      !!
46      !! ** Method  :   Set a gaussian anomaly of pressure and associated
47      !!                geostrophic velocities
48      !!----------------------------------------------------------------------
49      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(in   ) ::   pdept   ! depth of t-point               [m]
50      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(in   ) ::   ptmask  ! t-point ocean mask             [m]
51      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(  out) ::   pts     ! T & S fields      [Celsius ; g/kg]
52      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(  out) ::   pu      ! i-component of the velocity  [m/s]
53      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(  out) ::   pv      ! j-component of the velocity  [m/s]
54      REAL(wp), DIMENSION(jpi,jpj)         , INTENT(  out) ::   pssh    ! sea-surface height
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         ! sea level:
91         pssh(:,:) = 0.
92         ! temperature:
93         pts(:,:,:,jp_tem) = 10._wp
94         ! salinity: 
95         pts(:,:,:,jp_sal) = 35._wp
96         ! velocities:
97         pu(:,:,:) = 0.
98         pv(:,:,:) = 0.
99         
100      CASE(1)    ! geostrophic zonal jet from -zjety to +zjety
101
102         ! sea level:
103         SELECT CASE( nn_fcase )
104         CASE(0)    ! f = f0
105            ! sea level: ssh = - fuy / g
106            WHERE( ABS(gphit) <= zjety )
107               pssh(:,:) = - ff_t(:,:) * rn_uzonal * gphit(:,:) * 1.e3 / grav
108            ELSEWHERE
109               pssh(:,:) = - ff_t(:,:) * rn_uzonal * SIGN(zjety, gphit(:,:)) * 1.e3 / grav
110            END WHERE
111         CASE(1)    ! f = f0 + beta*y
112            ! sea level: ssh = - u / g * ( fy + 0.5 * beta * y^2 )
113            zbeta = 2._wp * omega * COS( rad * rn_ppgphi0 ) / ra
114            WHERE( ABS(gphit) <= zjety )
115               pssh(:,:) = - rn_uzonal / grav * ( zf0 * gphit(:,:) * 1.e3 + 0.5 * zbeta * gphit(:,:) * gphit(:,:) * 1.e6 )
116            ELSEWHERE
117               pssh(:,:) = - rn_uzonal / grav * ( zf0 * SIGN(zjety, gphit(:,:)) * 1.e3   &
118                  &                             + 0.5 * zbeta * zjety * zjety * 1.e6 )
119            END WHERE
120         END SELECT
121         ! temperature:
122         pts(:,:,:,jp_tem) = 10._wp
123         ! salinity: 
124         pts(:,:,jpk,jp_sal) = 0.
125         DO jk=1, jpkm1
126            WHERE( ABS(gphit) <= zjety )
127!!$            WHERE( ABS(gphit) <= zjety*0.5 .AND. ABS(glamt) <= zjety*0.5 ) ! for a square of salt
128               pts(:,:,jk,jp_sal) = 35.
129            ELSEWHERE
130               pts(:,:,jk,jp_sal) = 30.
131            END WHERE                   
132         END DO
133         ! velocities:
134         pu(:,:,:) = 0.
135         DO jk=1, jpkm1
136            WHERE( ABS(gphit) <= zjety ) pu(:,:,jk) = rn_uzonal
137         END DO
138         pv(:,:,:) = 0.
139         !                 
140      CASE(2)    ! geostrophic zonal current shear
141     
142         ! sea level:
143         SELECT CASE( nn_fcase )
144         CASE(0)    ! f = f0
145            ! sea level: ssh = - fuy / g
146            WHERE( ABS(gphit) <= zjety )
147               pssh(:,:) = - ff_t(:,:) * rn_uzonal * ABS(gphit(:,:)) * 1.e3 / grav
148            ELSEWHERE
149               pssh(:,:) = - ff_t(:,:) * rn_uzonal * zjety * 1.e3 / grav
150            END WHERE
151         CASE(1)    ! f = f0 + beta*y
152            ! sea level: ssh = - u / g * ( fy + 0.5 * beta * y^2 )
153            zbeta = 2._wp * omega * COS( rad * rn_ppgphi0 ) / ra
154            WHERE( ABS(gphit) <= zjety )
155               pssh(:,:) = - SIGN(rn_uzonal, gphit(:,:)) / grav   &
156                  &        * ( zf0 * gphit(:,:) * 1.e3 + 0.5 * zbeta * gphit(:,:) * gphit(:,:) * 1.e6 )
157            ELSEWHERE
158               pssh(:,:) = - SIGN(rn_uzonal, gphit(:,:)) / grav   &
159                  &        * ( zf0 * SIGN(zjety, gphit(:,:)) * 1.e3 + 0.5 * zbeta * zjety * zjety * 1.e6 )
160            END WHERE
161         END SELECT
162         ! temperature:
163         pts(:,:,:,jp_tem) = 10._wp
164         ! salinity: 
165         pts(:,:,:,jp_sal) = 30.
166         DO jk=1, jpkm1
167            WHERE( ABS(gphiv) <= zjety ) pts(:,:,jk,jp_sal) = 30. + SIGN(1.,gphiv(:,:))
168         END DO
169         ! velocities:
170         pu(:,:,:) = 0.
171         DO jk=1, jpkm1
172            WHERE( ABS(gphiv) <= zjety ) pu(:,:,jk) = SIGN(rn_uzonal,gphit(:,:))*SIGN(1.,rn_uzonal)
173            WHERE( ABS(gphiv) == 0.    ) pu(:,:,jk) = 0. 
174         END DO
175         pv(:,:,:) = 0.
176         !                 
177      CASE(3)    ! gaussian zonal currant
178
179         ! zonal current
180         DO jk=1, jpkm1
181            ! gphit and lambda are both in km
182            pu(:,:,jk) = rn_uzonal * EXP( - 0.5 * gphit(:,:)**2 / rn_lambda**2 )
183         END DO
184         
185         ! sea level:
186         pssh(:,1) = - ff_t(:,1) / grav * pu(:,1,1) * e2t(:,1)
187         DO jl=1, jpnj
188            DO_2D( 0, 0, 0, 0 )
189               pssh(ji,jj) = pssh(ji,jj-1) - ff_t(ji,jj) / grav * pu(ji,jj,1) * e2t(ji,jj)
190            END_2D
191            CALL lbc_lnk( 'usrdef_istate', pssh, 'T',  1. )
192         END DO
193         
194         ! temperature:
195         pts(:,:,:,jp_tem) = 10._wp
196         ! salinity: 
197         DO jk=1, jpkm1
198            pts(:,:,jk,jp_sal) = pssh(:,:)
199         END DO
200         ! velocities:
201         pv(:,:,:) = 0.
202         !           
203      CASE(4)    ! geostrophic zonal pulse
204   
205         DO_2D( 1, 1, 1, 1 )
206            IF ( ABS(glamt(ji,jj)) <= zjetx ) THEN
207               zdu = rn_uzonal
208            ELSEIF ( ABS(glamt(ji,jj)) <= zjetx + 100. ) THEN
209               zdu = rn_uzonal * ( ( zjetx-ABS(glamt(ji,jj)) )/100. + 1. )
210            ELSE
211               zdu = 0.
212            END IF
213            IF ( ABS(gphit(ji,jj)) <= zjety ) THEN
214               pssh(ji,jj) = - ff_t(ji,jj) * zdu * gphit(ji,jj) * 1.e3 / grav
215               pu(ji,jj,:) = zdu
216               pts(ji,jj,:,jp_sal) = zdu / rn_uzonal + 1.
217            ELSE
218               pssh(ji,jj) = - ff_t(ji,jj) * zdu * SIGN(zjety,gphit(ji,jj)) * 1.e3 / grav 
219               pu(ji,jj,:) = 0.
220               pts(ji,jj,:,jp_sal) = 1.
221            END IF
222         END_2D
223         
224         ! temperature:
225         pts(:,:,:,jp_tem) = 10._wp * ptmask(:,:,:)       
226         pv(:,:,:) = 0.
227         
228       CASE(5)    ! vortex
229                  !
230         zf0   = 2._wp * omega * SIN( rad * rn_ppgphi0 )
231         zumax = rn_vtxmax * SIGN(1._wp, zf0) ! Here Anticyclonic: set zumax=-1 for cyclonic
232         zlambda = SQRT(2._wp)*rn_lambda*1.e3       ! Horizontal scale in meters
233         zn2 = 3.e-3**2
234         zH = 0.5_wp * 5000._wp
235         !
236         zr_lambda2 = 1._wp / zlambda**2
237         zP0 = rho0 * zf0 * zumax * zlambda * SQRT(EXP(1._wp)/2._wp)
238         !
239         DO_2D( 1, 1, 1, 1 )
240            zx = glamt(ji,jj) * 1.e3
241            zy = gphit(ji,jj) * 1.e3
242            ! Surface pressure: P(x,y,z) = F(z) * Psurf(x,y)
243            zpsurf = zP0 * EXP(-(zx**2+zy**2)*zr_lambda2) - rho0 * ff_t(ji,jj) * rn_uzonal * zy
244            ! Sea level:
245            pssh(ji,jj) = 0.
246            DO jl=1,5
247               zdt = pssh(ji,jj)
248               zdzF = (1._wp - EXP(zdt-zH)) / (zH - 1._wp + EXP(-zH))   ! F'(z)
249               zrho1 = rho0 * (1._wp + zn2*zdt/grav) - zdzF * zpsurf / grav    ! -1/g Dz(P) = -1/g * F'(z) * Psurf(x,y)
250               pssh(ji,jj) = zpsurf / (zrho1*grav) * ptmask(ji,jj,1)   ! ssh = Psurf / (Rho*g)
251            END DO
252            ! temperature:
253            DO jk=1,jpk
254               zdt =  pdept(ji,jj,jk) 
255               zrho1 = rho0 * (1._wp + zn2*zdt/grav)
256               IF (zdt < zH) THEN
257                  zdzF = (1._wp-EXP(zdt-zH)) / (zH-1._wp + EXP(-zH))   ! F'(z)
258                  zrho1 = zrho1 - zdzF * zpsurf / grav    ! -1/g Dz(P) = -1/g * F'(z) * Psurf(x,y)
259               ENDIF
260               !               pts(ji,jj,jk,jp_tem) = (20._wp + (rho0-zrho1) / 0.28_wp) * ptmask(ji,jj,jk)
261               pts(ji,jj,jk,jp_tem) = (10._wp + (rho0-zrho1) / 0.28_wp) * ptmask(ji,jj,jk)
262            END DO
263         END_2D
264         !
265         ! salinity: 
266         pts(:,:,:,jp_sal) = 35._wp * ptmask(:,:,:) 
267         !
268         ! velocities:
269         za = 2._wp * zP0 / zlambda**2
270         DO_2D( 0, 0, 0, 0 )
271            zx = glamu(ji,jj) * 1.e3
272            zy = gphiu(ji,jj) * 1.e3
273            DO jk=1, jpk
274               zdu = 0.5_wp * (pdept(ji,jj,jk) + pdept(ji+1,jj,jk))
275               IF (zdu < zH) THEN
276                  zf = (zH-1._wp-zdu+EXP(zdu-zH)) / (zH-1._wp+EXP(-zH))
277                  zdyPs = - za * zy * EXP(-(zx**2+zy**2)*zr_lambda2) - rho0 * ff_t(ji,jj) * rn_uzonal
278                  pu(ji,jj,jk) = - zf / ( rho0 * ff_t(ji,jj) ) * zdyPs * ptmask(ji,jj,jk) * ptmask(ji+1,jj,jk)
279               ELSE
280                  pu(ji,jj,jk) = 0._wp
281               ENDIF
282            END DO
283         END_2D
284         !
285         DO_2D( 0, 0, 0, 0 )
286            zx = glamv(ji,jj) * 1.e3
287            zy = gphiv(ji,jj) * 1.e3
288            DO jk=1, jpk
289               zdv = 0.5_wp * (pdept(ji,jj,jk) + pdept(ji,jj+1,jk))
290               IF (zdv < zH) THEN
291                  zf = (zH-1._wp-zdv+EXP(zdv-zH)) / (zH-1._wp+EXP(-zH))
292                  zdxPs = - za * zx * EXP(-(zx**2+zy**2)*zr_lambda2)
293                  pv(ji,jj,jk) = zf / ( rho0 * ff_f(ji,jj) ) * zdxPs * ptmask(ji,jj,jk) * ptmask(ji,jj+1,jk)
294               ELSE
295                  pv(ji,jj,jk) = 0._wp
296               ENDIF
297            END DO
298         END_2D
299         !           
300      END SELECT
301     
302      IF (ln_sshnoise) THEN
303         CALL RANDOM_SEED()
304         CALL RANDOM_NUMBER(zrandom)
305         pssh(:,:) = pssh(:,:) + ( 0.1  * zrandom(:,:) - 0.05 )
306      END IF
307      CALL lbc_lnk( 'usrdef_istate', pssh, 'T',  1. )
308      CALL lbc_lnk( 'usrdef_istate', pts , 'T',  1. )
309      CALL lbc_lnk_multi( 'usrdef_istate', pu, 'U', -1., pv, 'V', -1. )
310
311   END SUBROUTINE usr_def_istate
312
313   !!======================================================================
314END MODULE usrdef_istate
Note: See TracBrowser for help on using the repository browser.