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.
eosbn2.F90 in NEMO/branches/2020/dev_14237_KERNEL-01_IMMERSE_SEAMOUNT/tests/SEAMOUNT/MY_SRC – NEMO

source: NEMO/branches/2020/dev_14237_KERNEL-01_IMMERSE_SEAMOUNT/tests/SEAMOUNT/MY_SRC/eosbn2.F90 @ 14380

Last change on this file since 14380 was 14380, checked in by ayoung, 3 years ago

Corrected implementation of exponential equation of state. See ticket #2480.

File size: 84.7 KB
Line 
1MODULE eosbn2
2   !!==============================================================================
3   !!                       ***  MODULE  eosbn2  ***
4   !! Equation Of Seawater : in situ density - Brunt-Vaisala frequency
5   !!==============================================================================
6   !! History :  OPA  ! 1989-03  (O. Marti)  Original code
7   !!            6.0  ! 1994-07  (G. Madec, M. Imbard)  add bn2
8   !!            6.0  ! 1994-08  (G. Madec)  Add Jackett & McDougall eos
9   !!            7.0  ! 1996-01  (G. Madec)  statement function for e3
10   !!            8.1  ! 1997-07  (G. Madec)  density instead of volumic mass
11   !!             -   ! 1999-02  (G. Madec, N. Grima) semi-implicit pressure gradient
12   !!            8.2  ! 2001-09  (M. Ben Jelloul)  bugfix on linear eos
13   !!   NEMO     1.0  ! 2002-10  (G. Madec)  add eos_init
14   !!             -   ! 2002-11  (G. Madec, A. Bozec)  partial step, eos_insitu_2d
15   !!             -   ! 2003-08  (G. Madec)  F90, free form
16   !!            3.0  ! 2006-08  (G. Madec)  add tfreez function (now eos_fzp function)
17   !!            3.3  ! 2010-05  (C. Ethe, G. Madec)  merge TRC-TRA
18   !!             -   ! 2010-10  (G. Nurser, G. Madec)  add alpha/beta used in ldfslp
19   !!            3.7  ! 2012-03  (F. Roquet, G. Madec)  add primitive of alpha and beta used in PE computation
20   !!             -   ! 2012-05  (F. Roquet)  add Vallis and original JM95 equation of state
21   !!             -   ! 2013-04  (F. Roquet, G. Madec)  add eos_rab, change bn2 computation and reorganize the module
22   !!             -   ! 2014-09  (F. Roquet)  add TEOS-10, S-EOS, and modify EOS-80
23   !!             -   ! 2015-06  (P.A. Bouttier) eos_fzp functions changed to subroutines for AGRIF
24   !!----------------------------------------------------------------------
25
26   !!----------------------------------------------------------------------
27   !!   eos           : generic interface of the equation of state
28   !!   eos_insitu    : Compute the in situ density
29   !!   eos_insitu_pot: Compute the insitu and surface referenced potential volumic mass
30   !!   eos_insitu_2d : Compute the in situ density for 2d fields
31   !!   bn2           : compute the Brunt-Vaisala frequency
32   !!   eos_pt_from_ct: compute the potential temperature from the Conservative Temperature
33   !!   eos_rab       : generic interface of in situ thermal/haline expansion ratio
34   !!   eos_rab_3d    : compute in situ thermal/haline expansion ratio
35   !!   eos_rab_2d    : compute in situ thermal/haline expansion ratio for 2d fields
36   !!   eos_fzp_2d    : freezing temperature for 2d fields
37   !!   eos_fzp_0d    : freezing temperature for scalar
38   !!   eos_init      : set eos parameters (namelist)
39   !!----------------------------------------------------------------------
40   USE dom_oce        ! ocean space and time domain
41   USE domutl, ONLY : is_tile
42   USE phycst         ! physical constants
43   USE stopar         ! Stochastic T/S fluctuations
44   USE stopts         ! Stochastic T/S fluctuations
45   !
46   USE in_out_manager ! I/O manager
47   USE lib_mpp        ! MPP library
48   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined)
49   USE prtctl         ! Print control
50   USE lbclnk         ! ocean lateral boundary conditions
51   USE timing         ! Timing
52   USE usrdef_nam
53
54   IMPLICIT NONE
55   PRIVATE
56
57   !                  !! * Interface
58   INTERFACE eos
59      MODULE PROCEDURE eos_insitu, eos_insitu_pot, eos_insitu_2d, eos_insitu_pot_2d
60   END INTERFACE
61   !
62   INTERFACE eos_rab
63      MODULE PROCEDURE rab_3d, rab_2d, rab_0d
64   END INTERFACE
65   !
66   INTERFACE eos_fzp
67      MODULE PROCEDURE eos_fzp_2d, eos_fzp_0d
68   END INTERFACE
69   !
70   PUBLIC   eos            ! called by step, istate, tranpc and zpsgrd modules
71   PUBLIC   bn2            ! called by step module
72   PUBLIC   eos_rab        ! called by ldfslp, zdfddm, trabbl
73   PUBLIC   eos_pt_from_ct ! called by sbcssm
74   PUBLIC   eos_fzp        ! called by traadv_cen2 and sbcice_... modules
75   PUBLIC   eos_pen        ! used for pe diagnostics in trdpen module
76   PUBLIC   eos_init       ! called by istate module
77
78   !                               !!** Namelist nameos **
79   LOGICAL , PUBLIC ::   ln_TEOS10
80   LOGICAL , PUBLIC ::   ln_EOS80
81   LOGICAL , PUBLIC ::   ln_SEOS
82   LOGICAL , PUBLIC ::   ln_EEOS   ! determine if exponential eos is used
83
84   ! Parameters
85   LOGICAL , PUBLIC    ::   l_useCT         ! =T in ln_TEOS10=T (i.e. use eos_pt_from_ct to compute sst_m), =F otherwise
86   INTEGER , PUBLIC    ::   neos            ! Identifier for equation of state used
87
88   INTEGER , PARAMETER ::   np_teos10 = -1  ! parameter for using TEOS10
89   INTEGER , PARAMETER ::   np_eos80  =  0  ! parameter for using EOS80
90   INTEGER , PARAMETER ::   np_seos   = 1   ! parameter for using Simplified Equation of state
91   INTEGER , PARAMETER ::   np_eeos   =  2  ! parameter for using exponential equation of state (SEAMOUNT)
92
93   !                               !!!  simplified eos coefficients (default value: Vallis 2006)
94   REAL(wp), PUBLIC ::   rn_a0      = 1.6550e-1_wp     ! thermal expansion coeff.
95   REAL(wp), PUBLIC ::   rn_b0      = 7.6554e-1_wp     ! saline  expansion coeff.
96   REAL(wp) ::   rn_lambda1 = 5.9520e-2_wp     ! cabbeling coeff. in T^2
97   REAL(wp) ::   rn_lambda2 = 5.4914e-4_wp     ! cabbeling coeff. in S^2
98   REAL(wp) ::   rn_mu1     = 1.4970e-4_wp     ! thermobaric coeff. in T
99   REAL(wp) ::   rn_mu2     = 1.1090e-5_wp     ! thermobaric coeff. in S
100   REAL(wp) ::   rn_nu      = 2.4341e-3_wp     ! cabbeling coeff. in theta*salt
101
102   ! TEOS10/EOS80 parameters
103   REAL(wp) ::   r1_S0, r1_T0, r1_Z0, rdeltaS
104
105   ! EOS parameters
106   REAL(wp) ::   EOS000 , EOS100 , EOS200 , EOS300 , EOS400 , EOS500 , EOS600
107   REAL(wp) ::   EOS010 , EOS110 , EOS210 , EOS310 , EOS410 , EOS510
108   REAL(wp) ::   EOS020 , EOS120 , EOS220 , EOS320 , EOS420
109   REAL(wp) ::   EOS030 , EOS130 , EOS230 , EOS330
110   REAL(wp) ::   EOS040 , EOS140 , EOS240
111   REAL(wp) ::   EOS050 , EOS150
112   REAL(wp) ::   EOS060
113   REAL(wp) ::   EOS001 , EOS101 , EOS201 , EOS301 , EOS401
114   REAL(wp) ::   EOS011 , EOS111 , EOS211 , EOS311
115   REAL(wp) ::   EOS021 , EOS121 , EOS221
116   REAL(wp) ::   EOS031 , EOS131
117   REAL(wp) ::   EOS041
118   REAL(wp) ::   EOS002 , EOS102 , EOS202
119   REAL(wp) ::   EOS012 , EOS112
120   REAL(wp) ::   EOS022
121   REAL(wp) ::   EOS003 , EOS103
122   REAL(wp) ::   EOS013
123
124   ! ALPHA parameters
125   REAL(wp) ::   ALP000 , ALP100 , ALP200 , ALP300 , ALP400 , ALP500
126   REAL(wp) ::   ALP010 , ALP110 , ALP210 , ALP310 , ALP410
127   REAL(wp) ::   ALP020 , ALP120 , ALP220 , ALP320
128   REAL(wp) ::   ALP030 , ALP130 , ALP230
129   REAL(wp) ::   ALP040 , ALP140
130   REAL(wp) ::   ALP050
131   REAL(wp) ::   ALP001 , ALP101 , ALP201 , ALP301
132   REAL(wp) ::   ALP011 , ALP111 , ALP211
133   REAL(wp) ::   ALP021 , ALP121
134   REAL(wp) ::   ALP031
135   REAL(wp) ::   ALP002 , ALP102
136   REAL(wp) ::   ALP012
137   REAL(wp) ::   ALP003
138
139   ! BETA parameters
140   REAL(wp) ::   BET000 , BET100 , BET200 , BET300 , BET400 , BET500
141   REAL(wp) ::   BET010 , BET110 , BET210 , BET310 , BET410
142   REAL(wp) ::   BET020 , BET120 , BET220 , BET320
143   REAL(wp) ::   BET030 , BET130 , BET230
144   REAL(wp) ::   BET040 , BET140
145   REAL(wp) ::   BET050
146   REAL(wp) ::   BET001 , BET101 , BET201 , BET301
147   REAL(wp) ::   BET011 , BET111 , BET211
148   REAL(wp) ::   BET021 , BET121
149   REAL(wp) ::   BET031
150   REAL(wp) ::   BET002 , BET102
151   REAL(wp) ::   BET012
152   REAL(wp) ::   BET003
153
154   ! PEN parameters
155   REAL(wp) ::   PEN000 , PEN100 , PEN200 , PEN300 , PEN400
156   REAL(wp) ::   PEN010 , PEN110 , PEN210 , PEN310
157   REAL(wp) ::   PEN020 , PEN120 , PEN220
158   REAL(wp) ::   PEN030 , PEN130
159   REAL(wp) ::   PEN040
160   REAL(wp) ::   PEN001 , PEN101 , PEN201
161   REAL(wp) ::   PEN011 , PEN111
162   REAL(wp) ::   PEN021
163   REAL(wp) ::   PEN002 , PEN102
164   REAL(wp) ::   PEN012
165
166   ! ALPHA_PEN parameters
167   REAL(wp) ::   APE000 , APE100 , APE200 , APE300
168   REAL(wp) ::   APE010 , APE110 , APE210
169   REAL(wp) ::   APE020 , APE120
170   REAL(wp) ::   APE030
171   REAL(wp) ::   APE001 , APE101
172   REAL(wp) ::   APE011
173   REAL(wp) ::   APE002
174
175   ! BETA_PEN parameters
176   REAL(wp) ::   BPE000 , BPE100 , BPE200 , BPE300
177   REAL(wp) ::   BPE010 , BPE110 , BPE210
178   REAL(wp) ::   BPE020 , BPE120
179   REAL(wp) ::   BPE030
180   REAL(wp) ::   BPE001 , BPE101
181   REAL(wp) ::   BPE011
182   REAL(wp) ::   BPE002
183
184   !! * Substitutions
185#  include "do_loop_substitute.h90"
186#  include "domzgr_substitute.h90"
187   !!----------------------------------------------------------------------
188   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
189   !! $Id: eosbn2.F90 14131 2020-12-08 12:09:14Z jchanut $
190   !! Software governed by the CeCILL license (see ./LICENSE)
191   !!----------------------------------------------------------------------
192CONTAINS
193
194   SUBROUTINE eos_insitu( pts, prd, pdep )
195      !!
196      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pts   ! 1 : potential temperature  [Celsius]
197      !                                                      ! 2 : salinity               [psu]
198      REAL(wp), DIMENSION(:,:,:)  , INTENT(  out) ::   prd   ! in situ density            [-]
199      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pdep  ! depth                      [m]
200      !!
201      CALL eos_insitu_t( pts, is_tile(pts), prd, is_tile(prd), pdep, is_tile(pdep) )
202   END SUBROUTINE eos_insitu
203
204   SUBROUTINE eos_insitu_t( pts, ktts, prd, ktrd, pdep, ktdep )
205      !!----------------------------------------------------------------------
206      !!                   ***  ROUTINE eos_insitu  ***
207      !!
208      !! ** Purpose :   Compute the in situ density (ratio rho/rho0) from
209      !!       potential temperature and salinity using an equation of state
210      !!       selected in the nameos namelist
211      !!
212      !! ** Method  :   prd(t,s,z) = ( rho(t,s,z) - rho0 ) / rho0
213      !!         with   prd    in situ density anomaly      no units
214      !!                t      TEOS10: CT or EOS80: PT      Celsius
215      !!                s      TEOS10: SA or EOS80: SP      TEOS10: g/kg or EOS80: psu
216      !!                z      depth                        meters
217      !!                rho    in situ density              kg/m^3
218      !!                rho0   reference density            kg/m^3
219      !!
220      !!     ln_teos10 : polynomial TEOS-10 equation of state is used for rho(t,s,z).
221      !!         Check value: rho = 1028.21993233072 kg/m^3 for z=3000 dbar, ct=3 Celsius, sa=35.5 g/kg
222      !!
223      !!     ln_eos80 : polynomial EOS-80 equation of state is used for rho(t,s,z).
224      !!         Check value: rho = 1028.35011066567 kg/m^3 for z=3000 dbar, pt=3 Celsius, sp=35.5 psu
225      !!
226      !!     ln_seos : simplified equation of state
227      !!              prd(t,s,z) = ( -a0*(1+lambda/2*(T-T0)+mu*z+nu*(S-S0))*(T-T0) + b0*(S-S0) ) / rho0
228      !!              linear case function of T only: rn_alpha<>0, other coefficients = 0
229      !!              linear eos function of T and S: rn_alpha and rn_beta<>0, other coefficients=0
230      !!              Vallis like equation: use default values of coefficients
231      !!     ln_eosano : exponential SEAMOUNT equation of state returning density anomaly only
232      !!              prd(t,s,z) = ( -a0*(T-T0) + b0*(S-S0) ) / rho0
233      !!              setup for SEAMOUNT exponential eos
234      !!
235      !! ** Action  :   compute prd , the in situ density (no units)
236      !!
237      !! References :   Roquet et al, Ocean Modelling, in preparation (2014)
238      !!                Vallis, Atmospheric and Oceanic Fluid Dynamics, 2006
239      !!                TEOS-10 Manual, 2010
240      !!----------------------------------------------------------------------
241      INTEGER                                 , INTENT(in   ) ::   ktts, ktrd, ktdep
242      REAL(wp), DIMENSION(A2D_T(ktts) ,JPK,JPTS), INTENT(in   ) ::   pts   ! 1 : potential temperature  [Celsius]
243      !                                                                  ! 2 : salinity               [psu]
244      REAL(wp), DIMENSION(A2D_T(ktrd) ,JPK     ), INTENT(  out) ::   prd   ! in situ density            [-]
245      REAL(wp), DIMENSION(A2D_T(ktdep),JPK     ), INTENT(in   ) ::   pdep  ! depth                      [m]
246      !
247      INTEGER  ::   ji, jj, jk                ! dummy loop indices
248      REAL(wp) ::   zt , zh , zs , ztm        ! local scalars
249      REAL(wp) ::   zn , zn0, zn1, zn2, zn3   !   -      -
250      !!----------------------------------------------------------------------
251      !
252      IF( ln_timing )   CALL timing_start('eos-insitu')
253      !
254      SELECT CASE( neos )
255      !
256      CASE( np_teos10, np_eos80 )                !==  polynomial TEOS-10 / EOS-80 ==!
257         !
258         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
259            !
260            zh  = pdep(ji,jj,jk) * r1_Z0                                  ! depth
261            zt  = pts (ji,jj,jk,jp_tem) * r1_T0                           ! temperature
262            zs  = SQRT( ABS( pts(ji,jj,jk,jp_sal) + rdeltaS ) * r1_S0 )   ! square root salinity
263            ztm = tmask(ji,jj,jk)                                         ! tmask
264            !
265            zn3 = EOS013*zt   &
266               &   + EOS103*zs+EOS003
267               !
268            zn2 = (EOS022*zt   &
269               &   + EOS112*zs+EOS012)*zt   &
270               &   + (EOS202*zs+EOS102)*zs+EOS002
271               !
272            zn1 = (((EOS041*zt   &
273               &   + EOS131*zs+EOS031)*zt   &
274               &   + (EOS221*zs+EOS121)*zs+EOS021)*zt   &
275               &   + ((EOS311*zs+EOS211)*zs+EOS111)*zs+EOS011)*zt   &
276               &   + (((EOS401*zs+EOS301)*zs+EOS201)*zs+EOS101)*zs+EOS001
277               !
278            zn0 = (((((EOS060*zt   &
279               &   + EOS150*zs+EOS050)*zt   &
280               &   + (EOS240*zs+EOS140)*zs+EOS040)*zt   &
281               &   + ((EOS330*zs+EOS230)*zs+EOS130)*zs+EOS030)*zt   &
282               &   + (((EOS420*zs+EOS320)*zs+EOS220)*zs+EOS120)*zs+EOS020)*zt   &
283               &   + ((((EOS510*zs+EOS410)*zs+EOS310)*zs+EOS210)*zs+EOS110)*zs+EOS010)*zt   &
284               &   + (((((EOS600*zs+EOS500)*zs+EOS400)*zs+EOS300)*zs+EOS200)*zs+EOS100)*zs+EOS000
285               !
286            zn  = ( ( zn3 * zh + zn2 ) * zh + zn1 ) * zh + zn0
287            !
288            prd(ji,jj,jk) = (  zn * r1_rho0 - 1._wp  ) * ztm  ! density anomaly (masked)
289            !
290         END_3D
291         !
292      CASE( np_seos )                !==  simplified EOS  ==!
293         !
294         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
295            zt  = pts  (ji,jj,jk,jp_tem) - 10._wp
296            zs  = pts  (ji,jj,jk,jp_sal) - 35._wp
297            zh  = pdep (ji,jj,jk)
298            ztm = tmask(ji,jj,jk)
299            !
300            zn =  - rn_a0 * ( 1._wp + 0.5_wp*rn_lambda1*zt + rn_mu1*zh ) * zt   &
301               &  + rn_b0 * ( 1._wp - 0.5_wp*rn_lambda2*zs - rn_mu2*zh ) * zs   &
302               &  - rn_nu * zt * zs
303               !
304            prd(ji,jj,jk) = zn * r1_rho0 * ztm                ! density anomaly (masked)
305         END_3D
306         !
307      CASE( np_eeos )                !==  exponential SEAMOUNT EOS  ==!
308         !
309         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
310            zt  = pts  (ji,jj,jk,jp_tem) + ( 28 - rn_drho * EXP(- pdep (ji,jj,jk) / 1000._wp) ) / rn_a0
311            ztm = tmask(ji,jj,jk)
312            !
313            zn = - rn_a0 * zt
314            !                                 
315            prd(ji,jj,jk) = zn * r1_rho0 * ztm                ! density anomaly (masked)
316         END_3D
317      END SELECT
318      !
319      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab3d_1=prd, clinfo1=' eos-insitu  : ', kdim=jpk )
320      !
321      IF( ln_timing )   CALL timing_stop('eos-insitu')
322      !
323   END SUBROUTINE eos_insitu_t
324
325
326   SUBROUTINE eos_insitu_pot( pts, prd, prhop, pdep )
327      !!
328      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pts    ! 1 : potential temperature  [Celsius]
329      !                                                       ! 2 : salinity               [psu]
330      REAL(wp), DIMENSION(:,:,:)  , INTENT(  out) ::   prd    ! in situ density            [-]
331      REAL(wp), DIMENSION(:,:,:)  , INTENT(  out) ::   prhop  ! potential density (surface referenced)
332      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pdep   ! depth                      [m]
333      !!
334      CALL eos_insitu_pot_t( pts, is_tile(pts), prd, is_tile(prd), prhop, is_tile(prhop), pdep, is_tile(pdep) )
335   END SUBROUTINE eos_insitu_pot
336
337
338   SUBROUTINE eos_insitu_pot_t( pts, ktts, prd, ktrd, prhop, ktrhop, pdep, ktdep )
339      !!----------------------------------------------------------------------
340      !!                  ***  ROUTINE eos_insitu_pot  ***
341      !!
342      !! ** Purpose :   Compute the in situ density (ratio rho/rho0) and the
343      !!      potential volumic mass (Kg/m3) from potential temperature and
344      !!      salinity fields using an equation of state selected in the
345      !!     namelist.
346      !!
347      !! ** Action  : - prd  , the in situ density (no units)
348      !!              - prhop, the potential volumic mass (Kg/m3)
349      !!
350      !!----------------------------------------------------------------------
351      INTEGER                                  , INTENT(in   ) ::   ktts, ktrd, ktrhop, ktdep
352      REAL(wp), DIMENSION(A2D_T(ktts)  ,JPK,JPTS), INTENT(in   ) ::   pts    ! 1 : potential temperature  [Celsius]
353      !                                                                    ! 2 : salinity               [psu]
354      REAL(wp), DIMENSION(A2D_T(ktrd)  ,JPK     ), INTENT(  out) ::   prd    ! in situ density            [-]
355      REAL(wp), DIMENSION(A2D_T(ktrhop),JPK     ), INTENT(  out) ::   prhop  ! potential density (surface referenced)
356      REAL(wp), DIMENSION(A2D_T(ktdep) ,JPK     ), INTENT(in   ) ::   pdep   ! depth                      [m]
357      !
358      INTEGER  ::   ji, jj, jk, jsmp             ! dummy loop indices
359      INTEGER  ::   jdof
360      REAL(wp) ::   zt , zh , zstemp, zs , ztm   ! local scalars
361      REAL(wp) ::   zn , zn0, zn1, zn2, zn3      !   -      -
362      REAL(wp), DIMENSION(:), ALLOCATABLE :: zn0_sto, zn_sto, zsign    ! local vectors
363      !!----------------------------------------------------------------------
364      !
365      IF( ln_timing )   CALL timing_start('eos-pot')
366      !
367      SELECT CASE ( neos )
368      !
369      CASE( np_teos10, np_eos80 )                !==  polynomial TEOS-10 / EOS-80 ==!
370         !
371         ! Stochastic equation of state
372         IF ( ln_sto_eos ) THEN
373            ALLOCATE(zn0_sto(1:2*nn_sto_eos))
374            ALLOCATE(zn_sto(1:2*nn_sto_eos))
375            ALLOCATE(zsign(1:2*nn_sto_eos))
376            DO jsmp = 1, 2*nn_sto_eos, 2
377              zsign(jsmp)   = 1._wp
378              zsign(jsmp+1) = -1._wp
379            END DO
380            !
381            DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
382               !
383               ! compute density (2*nn_sto_eos) times:
384               ! (1) for t+dt, s+ds (with the random TS fluctutation computed in sto_pts)
385               ! (2) for t-dt, s-ds (with the opposite fluctuation)
386               DO jsmp = 1, nn_sto_eos*2
387                  jdof   = (jsmp + 1) / 2
388                  zh     = pdep(ji,jj,jk) * r1_Z0                                  ! depth
389                  zt     = (pts (ji,jj,jk,jp_tem) + pts_ran(ji,jj,jk,jp_tem,jdof) * zsign(jsmp)) * r1_T0    ! temperature
390                  zstemp = pts  (ji,jj,jk,jp_sal) + pts_ran(ji,jj,jk,jp_sal,jdof) * zsign(jsmp)
391                  zs     = SQRT( ABS( zstemp + rdeltaS ) * r1_S0 )   ! square root salinity
392                  ztm    = tmask(ji,jj,jk)                                         ! tmask
393                  !
394                  zn3 = EOS013*zt   &
395                     &   + EOS103*zs+EOS003
396                     !
397                  zn2 = (EOS022*zt   &
398                     &   + EOS112*zs+EOS012)*zt   &
399                     &   + (EOS202*zs+EOS102)*zs+EOS002
400                     !
401                  zn1 = (((EOS041*zt   &
402                     &   + EOS131*zs+EOS031)*zt   &
403                     &   + (EOS221*zs+EOS121)*zs+EOS021)*zt   &
404                     &   + ((EOS311*zs+EOS211)*zs+EOS111)*zs+EOS011)*zt   &
405                     &   + (((EOS401*zs+EOS301)*zs+EOS201)*zs+EOS101)*zs+EOS001
406                     !
407                  zn0_sto(jsmp) = (((((EOS060*zt   &
408                     &   + EOS150*zs+EOS050)*zt   &
409                     &   + (EOS240*zs+EOS140)*zs+EOS040)*zt   &
410                     &   + ((EOS330*zs+EOS230)*zs+EOS130)*zs+EOS030)*zt   &
411                     &   + (((EOS420*zs+EOS320)*zs+EOS220)*zs+EOS120)*zs+EOS020)*zt   &
412                     &   + ((((EOS510*zs+EOS410)*zs+EOS310)*zs+EOS210)*zs+EOS110)*zs+EOS010)*zt   &
413                     &   + (((((EOS600*zs+EOS500)*zs+EOS400)*zs+EOS300)*zs+EOS200)*zs+EOS100)*zs+EOS000
414                     !
415                  zn_sto(jsmp)  = ( ( zn3 * zh + zn2 ) * zh + zn1 ) * zh + zn0_sto(jsmp)
416               END DO
417               !
418               ! compute stochastic density as the mean of the (2*nn_sto_eos) densities
419               prhop(ji,jj,jk) = 0._wp ; prd(ji,jj,jk) = 0._wp
420               DO jsmp = 1, nn_sto_eos*2
421                  prhop(ji,jj,jk) = prhop(ji,jj,jk) + zn0_sto(jsmp)                      ! potential density referenced at the surface
422                  !
423                  prd(ji,jj,jk) = prd(ji,jj,jk) + (  zn_sto(jsmp) * r1_rho0 - 1._wp  )   ! density anomaly (masked)
424               END DO
425               prhop(ji,jj,jk) = 0.5_wp * prhop(ji,jj,jk) * ztm / nn_sto_eos
426               prd  (ji,jj,jk) = 0.5_wp * prd  (ji,jj,jk) * ztm / nn_sto_eos
427            END_3D
428            DEALLOCATE(zn0_sto,zn_sto,zsign)
429         ! Non-stochastic equation of state
430         ELSE
431            DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
432               !
433               zh  = pdep(ji,jj,jk) * r1_Z0                                  ! depth
434               zt  = pts (ji,jj,jk,jp_tem) * r1_T0                           ! temperature
435               zs  = SQRT( ABS( pts(ji,jj,jk,jp_sal) + rdeltaS ) * r1_S0 )   ! square root salinity
436               ztm = tmask(ji,jj,jk)                                         ! tmask
437               !
438               zn3 = EOS013*zt   &
439                  &   + EOS103*zs+EOS003
440                  !
441               zn2 = (EOS022*zt   &
442                  &   + EOS112*zs+EOS012)*zt   &
443                  &   + (EOS202*zs+EOS102)*zs+EOS002
444                  !
445               zn1 = (((EOS041*zt   &
446                  &   + EOS131*zs+EOS031)*zt   &
447                  &   + (EOS221*zs+EOS121)*zs+EOS021)*zt   &
448                  &   + ((EOS311*zs+EOS211)*zs+EOS111)*zs+EOS011)*zt   &
449                  &   + (((EOS401*zs+EOS301)*zs+EOS201)*zs+EOS101)*zs+EOS001
450                  !
451               zn0 = (((((EOS060*zt   &
452                  &   + EOS150*zs+EOS050)*zt   &
453                  &   + (EOS240*zs+EOS140)*zs+EOS040)*zt   &
454                  &   + ((EOS330*zs+EOS230)*zs+EOS130)*zs+EOS030)*zt   &
455                  &   + (((EOS420*zs+EOS320)*zs+EOS220)*zs+EOS120)*zs+EOS020)*zt   &
456                  &   + ((((EOS510*zs+EOS410)*zs+EOS310)*zs+EOS210)*zs+EOS110)*zs+EOS010)*zt   &
457                  &   + (((((EOS600*zs+EOS500)*zs+EOS400)*zs+EOS300)*zs+EOS200)*zs+EOS100)*zs+EOS000
458                  !
459               zn  = ( ( zn3 * zh + zn2 ) * zh + zn1 ) * zh + zn0
460               !
461               prhop(ji,jj,jk) = zn0 * ztm                           ! potential density referenced at the surface
462               !
463               prd(ji,jj,jk) = (  zn * r1_rho0 - 1._wp  ) * ztm      ! density anomaly (masked)
464            END_3D
465         ENDIF
466
467      CASE( np_seos )                !==  simplified EOS  ==!
468         !
469         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
470            zt  = pts  (ji,jj,jk,jp_tem) - 10._wp
471            zs  = pts  (ji,jj,jk,jp_sal) - 35._wp
472            zh  = pdep (ji,jj,jk)
473            ztm = tmask(ji,jj,jk)
474            !                                                     ! potential density referenced at the surface
475            zn =  - rn_a0 * ( 1._wp + 0.5_wp*rn_lambda1*zt ) * zt   &
476               &  + rn_b0 * ( 1._wp - 0.5_wp*rn_lambda2*zs ) * zs   &
477               &  - rn_nu * zt * zs
478            prhop(ji,jj,jk) = ( rho0 + zn ) * ztm
479            !                                                     ! density anomaly (masked)
480            zn = zn - ( rn_a0 * rn_mu1 * zt + rn_b0 * rn_mu2 * zs ) * zh
481            prd(ji,jj,jk) = zn * r1_rho0 * ztm
482            !
483         END_3D
484         !
485      CASE( np_eeos )                !==  exponential SEAMOUNT EOS  ==!
486         !
487         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
488            zt  = pts  (ji,jj,jk,jp_tem) + ( 28 - rn_drho * EXP(- pdep (ji,jj,jk) / 1000._wp) ) / rn_a0
489            ztm = tmask(ji,jj,jk)
490            !
491            zn = - rn_a0 * zt
492            !                                                     ! potential density referenced at the surface
493            prhop(ji,jj,jk) = ( rho0 + zn ) * ztm
494            prd(ji,jj,jk) = zn * r1_rho0 * ztm                 ! density anomaly (masked)
495         END_3D
496         !
497      END SELECT
498      !
499      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab3d_1=prd, clinfo1=' eos-pot: ', &
500         &                                  tab3d_2=prhop, clinfo2=' pot : ', kdim=jpk )
501      !
502      IF( ln_timing )   CALL timing_stop('eos-pot')
503      !
504   END SUBROUTINE eos_insitu_pot_t
505
506
507   SUBROUTINE eos_insitu_2d( pts, pdep, prd )
508      !!
509      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pts   ! 1 : potential temperature  [Celsius]
510      !                                                    ! 2 : salinity               [psu]
511      REAL(wp), DIMENSION(:,:)  , INTENT(in   ) ::   pdep  ! depth                      [m]
512      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   prd   ! in situ density
513      !!
514      CALL eos_insitu_2d_t( pts, is_tile(pts), pdep, is_tile(pdep), prd, is_tile(prd) )
515   END SUBROUTINE eos_insitu_2d
516
517
518   SUBROUTINE eos_insitu_2d_t( pts, ktts, pdep, ktdep, prd, ktrd )
519      !!----------------------------------------------------------------------
520      !!                  ***  ROUTINE eos_insitu_2d  ***
521      !!
522      !! ** Purpose :   Compute the in situ density (ratio rho/rho0) from
523      !!      potential temperature and salinity using an equation of state
524      !!      selected in the nameos namelist. * 2D field case
525      !!
526      !! ** Action  : - prd , the in situ density (no units) (unmasked)
527      !!
528      !!----------------------------------------------------------------------
529      INTEGER                            , INTENT(in   ) ::   ktts, ktdep, ktrd
530      REAL(wp), DIMENSION(A2D_T(ktts),JPTS), INTENT(in   ) ::   pts   ! 1 : potential temperature  [Celsius]
531      !                                                             ! 2 : salinity               [psu]
532      REAL(wp), DIMENSION(A2D_T(ktdep)    ), INTENT(in   ) ::   pdep  ! depth                      [m]
533      REAL(wp), DIMENSION(A2D_T(ktrd)     ), INTENT(  out) ::   prd   ! in situ density
534      !
535      INTEGER  ::   ji, jj, jk                ! dummy loop indices
536      REAL(wp) ::   zt , zh , zs              ! local scalars
537      REAL(wp) ::   zn , zn0, zn1, zn2, zn3   !   -      -
538      !!----------------------------------------------------------------------
539      !
540      IF( ln_timing )   CALL timing_start('eos2d')
541      !
542      prd(:,:) = 0._wp
543      !
544      SELECT CASE( neos )
545      !
546      CASE( np_teos10, np_eos80 )                !==  polynomial TEOS-10 / EOS-80 ==!
547         !
548         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
549            !
550            zh  = pdep(ji,jj) * r1_Z0                                  ! depth
551            zt  = pts (ji,jj,jp_tem) * r1_T0                           ! temperature
552            zs  = SQRT( ABS( pts(ji,jj,jp_sal) + rdeltaS ) * r1_S0 )   ! square root salinity
553            !
554            zn3 = EOS013*zt   &
555               &   + EOS103*zs+EOS003
556               !
557            zn2 = (EOS022*zt   &
558               &   + EOS112*zs+EOS012)*zt   &
559               &   + (EOS202*zs+EOS102)*zs+EOS002
560               !
561            zn1 = (((EOS041*zt   &
562               &   + EOS131*zs+EOS031)*zt   &
563               &   + (EOS221*zs+EOS121)*zs+EOS021)*zt   &
564               &   + ((EOS311*zs+EOS211)*zs+EOS111)*zs+EOS011)*zt   &
565               &   + (((EOS401*zs+EOS301)*zs+EOS201)*zs+EOS101)*zs+EOS001
566               !
567            zn0 = (((((EOS060*zt   &
568               &   + EOS150*zs+EOS050)*zt   &
569               &   + (EOS240*zs+EOS140)*zs+EOS040)*zt   &
570               &   + ((EOS330*zs+EOS230)*zs+EOS130)*zs+EOS030)*zt   &
571               &   + (((EOS420*zs+EOS320)*zs+EOS220)*zs+EOS120)*zs+EOS020)*zt   &
572               &   + ((((EOS510*zs+EOS410)*zs+EOS310)*zs+EOS210)*zs+EOS110)*zs+EOS010)*zt   &
573               &   + (((((EOS600*zs+EOS500)*zs+EOS400)*zs+EOS300)*zs+EOS200)*zs+EOS100)*zs+EOS000
574               !
575            zn  = ( ( zn3 * zh + zn2 ) * zh + zn1 ) * zh + zn0
576            !
577            prd(ji,jj) = zn * r1_rho0 - 1._wp               ! unmasked in situ density anomaly
578            !
579         END_2D
580         !
581      CASE( np_seos )                !==  simplified EOS  ==!
582         !
583         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
584            !
585            zt    = pts  (ji,jj,jp_tem)  - 10._wp
586            zs    = pts  (ji,jj,jp_sal)  - 35._wp
587            zh    = pdep (ji,jj)                         ! depth at the partial step level
588            !
589            zn =  - rn_a0 * ( 1._wp + 0.5_wp*rn_lambda1*zt + rn_mu1*zh ) * zt   &
590               &  + rn_b0 * ( 1._wp - 0.5_wp*rn_lambda2*zs - rn_mu2*zh ) * zs   &
591               &  - rn_nu * zt * zs
592               !
593            prd(ji,jj) = zn * r1_rho0               ! unmasked in situ density anomaly
594            !
595         END_2D
596         !
597      CASE( np_eeos )                !==  exponential SEAMOUNT EOS  ==!
598         !
599         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
600            zt  = pts  (ji,jj,jp_tem) + ( 28 - rn_drho * EXP(- pdep (ji,jj) / 1000._wp) ) / rn_a0
601            !
602            zn = - rn_a0 * zt
603            !
604            prd(ji,jj) = zn * r1_rho0                 ! density anomaly (masked)
605         END_2D
606         !
607      END SELECT
608      !
609      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab2d_1=prd, clinfo1=' eos2d: ' )
610      !
611      IF( ln_timing )   CALL timing_stop('eos2d')
612      !
613   END SUBROUTINE eos_insitu_2d_t
614
615
616   SUBROUTINE eos_insitu_pot_2d( pts, prhop )
617      !!----------------------------------------------------------------------
618      !!                  ***  ROUTINE eos_insitu_pot  ***
619      !!
620      !! ** Purpose :   Compute the in situ density (ratio rho/rho0) and the
621      !!      potential volumic mass (Kg/m3) from potential temperature and
622      !!      salinity fields using an equation of state selected in the
623      !!     namelist.
624      !!
625      !! ** Action  :
626      !!              - prhop, the potential volumic mass (Kg/m3)
627      !!
628      !!----------------------------------------------------------------------
629      REAL(wp), DIMENSION(jpi,jpj,jpts), INTENT(in   ) ::   pts    ! 1 : potential temperature  [Celsius]
630      !                                                                ! 2 : salinity               [psu]
631      REAL(wp), DIMENSION(jpi,jpj     ), INTENT(  out) ::   prhop  ! potential density (surface referenced)
632      !
633      INTEGER  ::   ji, jj, jk, jsmp             ! dummy loop indices
634      INTEGER  ::   jdof
635      REAL(wp) ::   zt , zh , zstemp, zs , ztm   ! local scalars
636      REAL(wp) ::   zn , zn0, zn1, zn2, zn3      !   -      -
637      REAL(wp), DIMENSION(:), ALLOCATABLE :: zn0_sto, zn_sto, zsign    ! local vectors
638      !!----------------------------------------------------------------------
639      !
640      IF( ln_timing )   CALL timing_start('eos-pot')
641      !
642      SELECT CASE ( neos )
643      !
644      CASE( np_teos10, np_eos80 )                !==  polynomial TEOS-10 / EOS-80 ==!
645         !
646            DO_2D( 1, 1, 1, 1 )
647               !
648               zt  = pts (ji,jj,jp_tem) * r1_T0                           ! temperature
649               zs  = SQRT( ABS( pts(ji,jj,jp_sal) + rdeltaS ) * r1_S0 )   ! square root salinity
650               ztm = tmask(ji,jj,1)                                         ! tmask
651               !
652               zn0 = (((((EOS060*zt   &
653                  &   + EOS150*zs+EOS050)*zt   &
654                  &   + (EOS240*zs+EOS140)*zs+EOS040)*zt   &
655                  &   + ((EOS330*zs+EOS230)*zs+EOS130)*zs+EOS030)*zt   &
656                  &   + (((EOS420*zs+EOS320)*zs+EOS220)*zs+EOS120)*zs+EOS020)*zt   &
657                  &   + ((((EOS510*zs+EOS410)*zs+EOS310)*zs+EOS210)*zs+EOS110)*zs+EOS010)*zt   &
658                  &   + (((((EOS600*zs+EOS500)*zs+EOS400)*zs+EOS300)*zs+EOS200)*zs+EOS100)*zs+EOS000
659                  !
660               !
661               prhop(ji,jj) = zn0 * ztm                           ! potential density referenced at the surface
662               !
663            END_2D
664
665      CASE( np_seos )                !==  simplified EOS  ==!
666         !
667         DO_2D( 1, 1, 1, 1 )
668            zt  = pts  (ji,jj,jp_tem) - 10._wp
669            zs  = pts  (ji,jj,jp_sal) - 35._wp
670            ztm = tmask(ji,jj,1)
671            !                                                     ! potential density referenced at the surface
672            zn =  - rn_a0 * ( 1._wp + 0.5_wp*rn_lambda1*zt ) * zt   &
673               &  + rn_b0 * ( 1._wp - 0.5_wp*rn_lambda2*zs ) * zs   &
674               &  - rn_nu * zt * zs
675            prhop(ji,jj) = ( rho0 + zn ) * ztm
676            !
677         END_2D
678         !
679      CASE( np_eeos )                !==  exponential SEAMOUNT EOS  ==!
680         !
681         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
682            zt  = pts  (ji,jj,jp_tem) + ( 28 ) / rn_a0
683            ztm = tmask(ji,jj,1)
684            !
685            zn = - rn_a0 * zt
686            !
687            prhop(ji,jj) = ( rho0 + zn ) * ztm
688         END_2D
689      END SELECT
690      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab2d_1=prhop, clinfo1=' pot: ', kdim=1 )
691      !
692      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab2d_1=prhop, clinfo1=' eos-pot: ' )
693      !
694      IF( ln_timing )   CALL timing_stop('eos-pot')
695      !
696   END SUBROUTINE eos_insitu_pot_2d
697
698
699   SUBROUTINE rab_3d( pts, pab, Kmm )
700      !!
701      INTEGER                     , INTENT(in   ) ::   Kmm   ! time level index
702      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pts   ! pot. temperature & salinity
703      REAL(wp), DIMENSION(:,:,:,:), INTENT(  out) ::   pab   ! thermal/haline expansion ratio
704      !!
705      CALL rab_3d_t( pts, is_tile(pts), pab, is_tile(pab), Kmm )
706   END SUBROUTINE rab_3d
707
708
709   SUBROUTINE rab_3d_t( pts, ktts, pab, ktab, Kmm )
710      !!----------------------------------------------------------------------
711      !!                 ***  ROUTINE rab_3d  ***
712      !!
713      !! ** Purpose :   Calculates thermal/haline expansion ratio at T-points
714      !!
715      !! ** Method  :   calculates alpha / beta at T-points
716      !!
717      !! ** Action  : - pab     : thermal/haline expansion ratio at T-points
718      !!----------------------------------------------------------------------
719      INTEGER                                , INTENT(in   ) ::   Kmm   ! time level index
720      INTEGER                                , INTENT(in   ) ::   ktts, ktab
721      REAL(wp), DIMENSION(A2D_T(ktts),JPK,JPTS), INTENT(in   ) ::   pts   ! pot. temperature & salinity
722      REAL(wp), DIMENSION(A2D_T(ktab),JPK,JPTS), INTENT(  out) ::   pab   ! thermal/haline expansion ratio
723      !
724      INTEGER  ::   ji, jj, jk                ! dummy loop indices
725      REAL(wp) ::   zt , zh , zs , ztm        ! local scalars
726      REAL(wp) ::   zn , zn0, zn1, zn2, zn3   !   -      -
727      !!----------------------------------------------------------------------
728      !
729      IF( ln_timing )   CALL timing_start('rab_3d')
730      !
731      SELECT CASE ( neos )
732      !
733      CASE( np_teos10, np_eos80 )                !==  polynomial TEOS-10 / EOS-80 ==!
734         !
735         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
736            !
737            zh  = gdept(ji,jj,jk,Kmm) * r1_Z0                                ! depth
738            zt  = pts (ji,jj,jk,jp_tem) * r1_T0                           ! temperature
739            zs  = SQRT( ABS( pts(ji,jj,jk,jp_sal) + rdeltaS ) * r1_S0 )   ! square root salinity
740            ztm = tmask(ji,jj,jk)                                         ! tmask
741            !
742            ! alpha
743            zn3 = ALP003
744            !
745            zn2 = ALP012*zt + ALP102*zs+ALP002
746            !
747            zn1 = ((ALP031*zt   &
748               &   + ALP121*zs+ALP021)*zt   &
749               &   + (ALP211*zs+ALP111)*zs+ALP011)*zt   &
750               &   + ((ALP301*zs+ALP201)*zs+ALP101)*zs+ALP001
751               !
752            zn0 = ((((ALP050*zt   &
753               &   + ALP140*zs+ALP040)*zt   &
754               &   + (ALP230*zs+ALP130)*zs+ALP030)*zt   &
755               &   + ((ALP320*zs+ALP220)*zs+ALP120)*zs+ALP020)*zt   &
756               &   + (((ALP410*zs+ALP310)*zs+ALP210)*zs+ALP110)*zs+ALP010)*zt   &
757               &   + ((((ALP500*zs+ALP400)*zs+ALP300)*zs+ALP200)*zs+ALP100)*zs+ALP000
758               !
759            zn  = ( ( zn3 * zh + zn2 ) * zh + zn1 ) * zh + zn0
760            !
761            pab(ji,jj,jk,jp_tem) = zn * r1_rho0 * ztm
762            !
763            ! beta
764            zn3 = BET003
765            !
766            zn2 = BET012*zt + BET102*zs+BET002
767            !
768            zn1 = ((BET031*zt   &
769               &   + BET121*zs+BET021)*zt   &
770               &   + (BET211*zs+BET111)*zs+BET011)*zt   &
771               &   + ((BET301*zs+BET201)*zs+BET101)*zs+BET001
772               !
773            zn0 = ((((BET050*zt   &
774               &   + BET140*zs+BET040)*zt   &
775               &   + (BET230*zs+BET130)*zs+BET030)*zt   &
776               &   + ((BET320*zs+BET220)*zs+BET120)*zs+BET020)*zt   &
777               &   + (((BET410*zs+BET310)*zs+BET210)*zs+BET110)*zs+BET010)*zt   &
778               &   + ((((BET500*zs+BET400)*zs+BET300)*zs+BET200)*zs+BET100)*zs+BET000
779               !
780            zn  = ( ( zn3 * zh + zn2 ) * zh + zn1 ) * zh + zn0
781            !
782            pab(ji,jj,jk,jp_sal) = zn / zs * r1_rho0 * ztm
783            !
784         END_3D
785         !
786      CASE( np_seos )                  !==  simplified EOS  ==!
787         !
788         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
789            zt  = pts (ji,jj,jk,jp_tem) - 10._wp   ! pot. temperature anomaly (t-T0)
790            zs  = pts (ji,jj,jk,jp_sal) - 35._wp   ! abs. salinity anomaly (s-S0)
791            zh  = gdept(ji,jj,jk,Kmm)                ! depth in meters at t-point
792            ztm = tmask(ji,jj,jk)                  ! land/sea bottom mask = surf. mask
793            !
794            zn  = rn_a0 * ( 1._wp + rn_lambda1*zt + rn_mu1*zh ) + rn_nu*zs
795            pab(ji,jj,jk,jp_tem) = zn * r1_rho0 * ztm   ! alpha
796            !
797            zn  = rn_b0 * ( 1._wp - rn_lambda2*zs - rn_mu2*zh ) - rn_nu*zt
798            pab(ji,jj,jk,jp_sal) = zn * r1_rho0 * ztm   ! beta
799            !
800         END_3D
801         !
802      CASE( np_eeos )                !==  exponential SEAMOUNT EOS  ==!
803         !
804         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
805            zt  = pts  (ji,jj,jk,jp_tem) + ( 28 - rn_drho * EXP(- gdept (ji,jj,jk,Kmm) / 1000._wp) ) / rn_a0
806            ztm = tmask(ji,jj,jk)
807            !
808            zn = rn_a0
809            !                                 
810            pab(ji,jj,jk,jp_tem) = zn * r1_rho0 * ztm   ! alpha
811            pab(ji,jj,jk,jp_sal) = 0
812            !
813         END_3D
814      CASE DEFAULT
815         WRITE(ctmp1,*) '          bad flag value for neos = ', neos
816         CALL ctl_stop( 'rab_3d:', ctmp1 )
817         !
818      END SELECT
819      !
820      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab3d_1=pab(:,:,:,jp_tem), clinfo1=' rab_3d_t: ', &
821         &                                  tab3d_2=pab(:,:,:,jp_sal), clinfo2=' rab_3d_s : ', kdim=jpk )
822      !
823      IF( ln_timing )   CALL timing_stop('rab_3d')
824      !
825   END SUBROUTINE rab_3d_t
826
827
828   SUBROUTINE rab_2d( pts, pdep, pab, Kmm )
829      !!
830      INTEGER                   , INTENT(in   ) ::   Kmm   ! time level index
831      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pts    ! pot. temperature & salinity
832      REAL(wp), DIMENSION(:,:)  , INTENT(in   ) ::   pdep   ! depth                  [m]
833      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   pab    ! thermal/haline expansion ratio
834      !!
835      CALL rab_2d_t(pts, is_tile(pts), pdep, is_tile(pdep), pab, is_tile(pab), Kmm)
836   END SUBROUTINE rab_2d
837
838
839   SUBROUTINE rab_2d_t( pts, ktts, pdep, ktdep, pab, ktab, Kmm )
840      !!----------------------------------------------------------------------
841      !!                 ***  ROUTINE rab_2d  ***
842      !!
843      !! ** Purpose :   Calculates thermal/haline expansion ratio for a 2d field (unmasked)
844      !!
845      !! ** Action  : - pab     : thermal/haline expansion ratio at T-points
846      !!----------------------------------------------------------------------
847      INTEGER                            , INTENT(in   ) ::   Kmm   ! time level index
848      INTEGER                            , INTENT(in   ) ::   ktts, ktdep, ktab
849      REAL(wp), DIMENSION(A2D_T(ktts),JPTS), INTENT(in   ) ::   pts    ! pot. temperature & salinity
850      REAL(wp), DIMENSION(A2D_T(ktdep)    ), INTENT(in   ) ::   pdep   ! depth                  [m]
851      REAL(wp), DIMENSION(A2D_T(ktab),JPTS), INTENT(  out) ::   pab    ! thermal/haline expansion ratio
852      !
853      INTEGER  ::   ji, jj, jk                ! dummy loop indices
854      REAL(wp) ::   zt , zh , zs              ! local scalars
855      REAL(wp) ::   zn , zn0, zn1, zn2, zn3   !   -      -
856      !!----------------------------------------------------------------------
857      !
858      IF( ln_timing )   CALL timing_start('rab_2d')
859      !
860      pab(:,:,:) = 0._wp
861      !
862      SELECT CASE ( neos )
863      !
864      CASE( np_teos10, np_eos80 )                !==  polynomial TEOS-10 / EOS-80 ==!
865         !
866         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
867            !
868            zh  = pdep(ji,jj) * r1_Z0                                  ! depth
869            zt  = pts (ji,jj,jp_tem) * r1_T0                           ! temperature
870            zs  = SQRT( ABS( pts(ji,jj,jp_sal) + rdeltaS ) * r1_S0 )   ! square root salinity
871            !
872            ! alpha
873            zn3 = ALP003
874            !
875            zn2 = ALP012*zt + ALP102*zs+ALP002
876            !
877            zn1 = ((ALP031*zt   &
878               &   + ALP121*zs+ALP021)*zt   &
879               &   + (ALP211*zs+ALP111)*zs+ALP011)*zt   &
880               &   + ((ALP301*zs+ALP201)*zs+ALP101)*zs+ALP001
881               !
882            zn0 = ((((ALP050*zt   &
883               &   + ALP140*zs+ALP040)*zt   &
884               &   + (ALP230*zs+ALP130)*zs+ALP030)*zt   &
885               &   + ((ALP320*zs+ALP220)*zs+ALP120)*zs+ALP020)*zt   &
886               &   + (((ALP410*zs+ALP310)*zs+ALP210)*zs+ALP110)*zs+ALP010)*zt   &
887               &   + ((((ALP500*zs+ALP400)*zs+ALP300)*zs+ALP200)*zs+ALP100)*zs+ALP000
888               !
889            zn  = ( ( zn3 * zh + zn2 ) * zh + zn1 ) * zh + zn0
890            !
891            pab(ji,jj,jp_tem) = zn * r1_rho0
892            !
893            ! beta
894            zn3 = BET003
895            !
896            zn2 = BET012*zt + BET102*zs+BET002
897            !
898            zn1 = ((BET031*zt   &
899               &   + BET121*zs+BET021)*zt   &
900               &   + (BET211*zs+BET111)*zs+BET011)*zt   &
901               &   + ((BET301*zs+BET201)*zs+BET101)*zs+BET001
902               !
903            zn0 = ((((BET050*zt   &
904               &   + BET140*zs+BET040)*zt   &
905               &   + (BET230*zs+BET130)*zs+BET030)*zt   &
906               &   + ((BET320*zs+BET220)*zs+BET120)*zs+BET020)*zt   &
907               &   + (((BET410*zs+BET310)*zs+BET210)*zs+BET110)*zs+BET010)*zt   &
908               &   + ((((BET500*zs+BET400)*zs+BET300)*zs+BET200)*zs+BET100)*zs+BET000
909               !
910            zn  = ( ( zn3 * zh + zn2 ) * zh + zn1 ) * zh + zn0
911            !
912            pab(ji,jj,jp_sal) = zn / zs * r1_rho0
913            !
914            !
915         END_2D
916         !
917      CASE( np_seos )                  !==  simplified EOS  ==!
918         !
919         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
920            !
921            zt    = pts  (ji,jj,jp_tem) - 10._wp   ! pot. temperature anomaly (t-T0)
922            zs    = pts  (ji,jj,jp_sal) - 35._wp   ! abs. salinity anomaly (s-S0)
923            zh    = pdep (ji,jj)                   ! depth at the partial step level
924            !
925            zn  = rn_a0 * ( 1._wp + rn_lambda1*zt + rn_mu1*zh ) + rn_nu*zs
926            pab(ji,jj,jp_tem) = zn * r1_rho0   ! alpha
927            !
928            zn  = rn_b0 * ( 1._wp - rn_lambda2*zs - rn_mu2*zh ) - rn_nu*zt
929            pab(ji,jj,jp_sal) = zn * r1_rho0   ! beta
930            !
931         END_2D
932         !
933      CASE( np_eeos )                !==  exponential SEAMOUNT EOS  ==!
934         !
935         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
936            zt  = pts  (ji,jj,jp_tem) + ( 28 - rn_drho * EXP(- pdep (ji,jj) / 1000._wp) ) / rn_a0
937            !
938            zn = rn_a0
939            !                                 
940            pab(ji,jj,jp_tem) = zn * r1_rho0   ! alpha
941            pab(ji,jj,jp_sal) = 0.0_wp
942         END_2D
943      CASE DEFAULT
944         WRITE(ctmp1,*) '          bad flag value for neos = ', neos
945         CALL ctl_stop( 'rab_2d:', ctmp1 )
946         !
947      END SELECT
948      !
949      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab2d_1=pab(:,:,jp_tem), clinfo1=' rab_2d_t: ', &
950         &                                  tab2d_2=pab(:,:,jp_sal), clinfo2=' rab_2d_s : ' )
951      !
952      IF( ln_timing )   CALL timing_stop('rab_2d')
953      !
954   END SUBROUTINE rab_2d_t
955
956
957   SUBROUTINE rab_0d( pts, pdep, pab, Kmm )
958      !!----------------------------------------------------------------------
959      !!                 ***  ROUTINE rab_0d  ***
960      !!
961      !! ** Purpose :   Calculates thermal/haline expansion ratio for a 2d field (unmasked)
962      !!
963      !! ** Action  : - pab     : thermal/haline expansion ratio at T-points
964      !!----------------------------------------------------------------------
965      INTEGER                              , INTENT(in   ) ::   Kmm   ! time level index
966      REAL(wp), DIMENSION(jpts)    , INTENT(in   ) ::   pts    ! pot. temperature & salinity
967      REAL(wp),                      INTENT(in   ) ::   pdep   ! depth                  [m]
968      REAL(wp), DIMENSION(jpts)    , INTENT(  out) ::   pab    ! thermal/haline expansion ratio
969      !
970      REAL(wp) ::   zt , zh , zs              ! local scalars
971      REAL(wp) ::   zn , zn0, zn1, zn2, zn3   !   -      -
972      !!----------------------------------------------------------------------
973      !
974      IF( ln_timing )   CALL timing_start('rab_0d')
975      !
976      pab(:) = 0._wp
977      !
978      SELECT CASE ( neos )
979      !
980      CASE( np_teos10, np_eos80 )      !==  polynomial TEOS-10 / EOS-80 ==!
981         !
982         !
983         zh  = pdep * r1_Z0                                  ! depth
984         zt  = pts (jp_tem) * r1_T0                           ! temperature
985         zs  = SQRT( ABS( pts(jp_sal) + rdeltaS ) * r1_S0 )   ! square root salinity
986         !
987         ! alpha
988         zn3 = ALP003
989         !
990         zn2 = ALP012*zt + ALP102*zs+ALP002
991         !
992         zn1 = ((ALP031*zt   &
993            &   + ALP121*zs+ALP021)*zt   &
994            &   + (ALP211*zs+ALP111)*zs+ALP011)*zt   &
995            &   + ((ALP301*zs+ALP201)*zs+ALP101)*zs+ALP001
996            !
997         zn0 = ((((ALP050*zt   &
998            &   + ALP140*zs+ALP040)*zt   &
999            &   + (ALP230*zs+ALP130)*zs+ALP030)*zt   &
1000            &   + ((ALP320*zs+ALP220)*zs+ALP120)*zs+ALP020)*zt   &
1001            &   + (((ALP410*zs+ALP310)*zs+ALP210)*zs+ALP110)*zs+ALP010)*zt   &
1002            &   + ((((ALP500*zs+ALP400)*zs+ALP300)*zs+ALP200)*zs+ALP100)*zs+ALP000
1003            !
1004         zn  = ( ( zn3 * zh + zn2 ) * zh + zn1 ) * zh + zn0
1005         !
1006         pab(jp_tem) = zn * r1_rho0
1007         !
1008         ! beta
1009         zn3 = BET003
1010         !
1011         zn2 = BET012*zt + BET102*zs+BET002
1012         !
1013         zn1 = ((BET031*zt   &
1014            &   + BET121*zs+BET021)*zt   &
1015            &   + (BET211*zs+BET111)*zs+BET011)*zt   &
1016            &   + ((BET301*zs+BET201)*zs+BET101)*zs+BET001
1017            !
1018         zn0 = ((((BET050*zt   &
1019            &   + BET140*zs+BET040)*zt   &
1020            &   + (BET230*zs+BET130)*zs+BET030)*zt   &
1021            &   + ((BET320*zs+BET220)*zs+BET120)*zs+BET020)*zt   &
1022            &   + (((BET410*zs+BET310)*zs+BET210)*zs+BET110)*zs+BET010)*zt   &
1023            &   + ((((BET500*zs+BET400)*zs+BET300)*zs+BET200)*zs+BET100)*zs+BET000
1024            !
1025         zn  = ( ( zn3 * zh + zn2 ) * zh + zn1 ) * zh + zn0
1026         !
1027         pab(jp_sal) = zn / zs * r1_rho0
1028         !
1029         !
1030         !
1031      CASE( np_seos )                  !==  simplified EOS  ==!
1032         !
1033         zt    = pts(jp_tem) - 10._wp   ! pot. temperature anomaly (t-T0)
1034         zs    = pts(jp_sal) - 35._wp   ! abs. salinity anomaly (s-S0)
1035         zh    = pdep                   ! depth at the partial step level
1036         !
1037         zn  = rn_a0 * ( 1._wp + rn_lambda1*zt + rn_mu1*zh ) + rn_nu*zs
1038         pab(jp_tem) = zn * r1_rho0   ! alpha
1039         !
1040         zn  = rn_b0 * ( 1._wp - rn_lambda2*zs - rn_mu2*zh ) - rn_nu*zt
1041         pab(jp_sal) = zn * r1_rho0   ! beta
1042         !
1043      CASE( np_eeos )                !==  exponential SEAMOUNT EOS  ==!
1044         !
1045         zt  = pts  (jp_tem) + ( 28 - rn_drho * EXP(- pdep / 1000._wp) ) / rn_a0
1046         !
1047         zn =   rn_a0
1048         !                                 
1049         pab(jp_tem) = zn * r1_rho0    ! alpha
1050         pab(jp_sal) = zn * r1_rho0    ! beta
1051!
1052      CASE DEFAULT
1053         WRITE(ctmp1,*) '          bad flag value for neos = ', neos
1054         CALL ctl_stop( 'rab_0d:', ctmp1 )
1055         !
1056      END SELECT
1057      !
1058      IF( ln_timing )   CALL timing_stop('rab_0d')
1059      !
1060   END SUBROUTINE rab_0d
1061
1062
1063   SUBROUTINE bn2( pts, pab, pn2, Kmm )
1064      !!
1065      INTEGER                              , INTENT(in   ) ::  Kmm   ! time level index
1066      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in   ) ::  pts   ! pot. temperature and salinity   [Celsius,psu]
1067      REAL(wp), DIMENSION(:,:,:,:)         , INTENT(in   ) ::  pab   ! thermal/haline expansion coef.  [Celsius-1,psu-1]
1068      REAL(wp), DIMENSION(:,:,:)           , INTENT(  out) ::  pn2   ! Brunt-Vaisala frequency squared [1/s^2]
1069      !!
1070      CALL bn2_t( pts, pab, is_tile(pab), pn2, is_tile(pn2), Kmm )
1071   END SUBROUTINE bn2
1072
1073
1074   SUBROUTINE bn2_t( pts, pab, ktab, pn2, ktn2, Kmm )
1075      !!----------------------------------------------------------------------
1076      !!                  ***  ROUTINE bn2  ***
1077      !!
1078      !! ** Purpose :   Compute the local Brunt-Vaisala frequency at the
1079      !!                time-step of the input arguments
1080      !!
1081      !! ** Method  :   pn2 = grav * (alpha dk[T] + beta dk[S] ) / e3w
1082      !!      where alpha and beta are given in pab, and computed on T-points.
1083      !!      N.B. N^2 is set one for all to zero at jk=1 in istate module.
1084      !!
1085      !! ** Action  :   pn2 : square of the brunt-vaisala frequency at w-point
1086      !!
1087      !!----------------------------------------------------------------------
1088      INTEGER                                , INTENT(in   ) ::  Kmm   ! time level index
1089      INTEGER                                , INTENT(in   ) ::  ktab, ktn2
1090      REAL(wp), DIMENSION(jpi,jpj,  jpk,jpts), INTENT(in   ) ::  pts   ! pot. temperature and salinity   [Celsius,psu]
1091      REAL(wp), DIMENSION(A2D_T(ktab),JPK,JPTS), INTENT(in   ) ::  pab   ! thermal/haline expansion coef.  [Celsius-1,psu-1]
1092      REAL(wp), DIMENSION(A2D_T(ktn2),JPK     ), INTENT(  out) ::  pn2   ! Brunt-Vaisala frequency squared [1/s^2]
1093      !
1094      INTEGER  ::   ji, jj, jk      ! dummy loop indices
1095      REAL(wp) ::   zaw, zbw, zrw   ! local scalars
1096      !!----------------------------------------------------------------------
1097      !
1098      IF( ln_timing )   CALL timing_start('bn2')
1099      !
1100      DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 2, jpkm1 )      ! interior points only (2=< jk =< jpkm1 ); surface and bottom value set to zero one for all in istate.F90
1101         zrw =   ( gdepw(ji,jj,jk  ,Kmm) - gdept(ji,jj,jk,Kmm) )   &
1102            &  / ( gdept(ji,jj,jk-1,Kmm) - gdept(ji,jj,jk,Kmm) )
1103            !
1104         zaw = pab(ji,jj,jk,jp_tem) * (1. - zrw) + pab(ji,jj,jk-1,jp_tem) * zrw
1105         zbw = pab(ji,jj,jk,jp_sal) * (1. - zrw) + pab(ji,jj,jk-1,jp_sal) * zrw
1106         !
1107         pn2(ji,jj,jk) = grav * (  zaw * ( pts(ji,jj,jk-1,jp_tem) - pts(ji,jj,jk,jp_tem) )     &
1108            &                    - zbw * ( pts(ji,jj,jk-1,jp_sal) - pts(ji,jj,jk,jp_sal) )  )  &
1109            &            / e3w(ji,jj,jk,Kmm) * wmask(ji,jj,jk)
1110      END_3D
1111      !
1112      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab3d_1=pn2, clinfo1=' bn2  : ', kdim=jpk )
1113      !
1114      IF( ln_timing )   CALL timing_stop('bn2')
1115      !
1116   END SUBROUTINE bn2_t
1117
1118
1119   FUNCTION eos_pt_from_ct( ctmp, psal ) RESULT( ptmp )
1120      !!----------------------------------------------------------------------
1121      !!                 ***  ROUTINE eos_pt_from_ct  ***
1122      !!
1123      !! ** Purpose :   Compute pot.temp. from cons. temp. [Celsius]
1124      !!
1125      !! ** Method  :   rational approximation (5/3th order) of TEOS-10 algorithm
1126      !!       checkvalue: pt=20.02391895 Celsius for sa=35.7g/kg, ct=20degC
1127      !!
1128      !! Reference  :   TEOS-10, UNESCO
1129      !!                Rational approximation to TEOS10 algorithm (rms error on WOA13 values: 4.0e-5 degC)
1130      !!----------------------------------------------------------------------
1131      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ) ::   ctmp   ! Cons. Temp   [Celsius]
1132      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ) ::   psal   ! salinity     [psu]
1133      ! Leave result array automatic rather than making explicitly allocated
1134      REAL(wp), DIMENSION(jpi,jpj) ::   ptmp   ! potential temperature [Celsius]
1135      !
1136      INTEGER  ::   ji, jj               ! dummy loop indices
1137      REAL(wp) ::   zt , zs , ztm        ! local scalars
1138      REAL(wp) ::   zn , zd              ! local scalars
1139      REAL(wp) ::   zdeltaS , z1_S0 , z1_T0
1140      !!----------------------------------------------------------------------
1141      !
1142      IF( ln_timing )   CALL timing_start('eos_pt_from_ct')
1143      !
1144      zdeltaS = 5._wp
1145      z1_S0   = 0.875_wp/35.16504_wp
1146      z1_T0   = 1._wp/40._wp
1147      !
1148      DO_2D( 1, 1, 1, 1 )
1149         !
1150         zt  = ctmp   (ji,jj) * z1_T0
1151         zs  = SQRT( ABS( psal(ji,jj) + zdeltaS ) * r1_S0 )
1152         ztm = tmask(ji,jj,1)
1153         !
1154         zn = ((((-2.1385727895e-01_wp*zt   &
1155            &   - 2.7674419971e-01_wp*zs+1.0728094330_wp)*zt   &
1156            &   + (2.6366564313_wp*zs+3.3546960647_wp)*zs-7.8012209473_wp)*zt   &
1157            &   + ((1.8835586562_wp*zs+7.3949191679_wp)*zs-3.3937395875_wp)*zs-5.6414948432_wp)*zt   &
1158            &   + (((3.5737370589_wp*zs-1.5512427389e+01_wp)*zs+2.4625741105e+01_wp)*zs   &
1159            &      +1.9912291000e+01_wp)*zs-3.2191146312e+01_wp)*zt   &
1160            &   + ((((5.7153204649e-01_wp*zs-3.0943149543_wp)*zs+9.3052495181_wp)*zs   &
1161            &      -9.4528934807_wp)*zs+3.1066408996_wp)*zs-4.3504021262e-01_wp
1162            !
1163         zd = (2.0035003456_wp*zt   &
1164            &   -3.4570358592e-01_wp*zs+5.6471810638_wp)*zt   &
1165            &   + (1.5393993508_wp*zs-6.9394762624_wp)*zs+1.2750522650e+01_wp
1166            !
1167         ptmp(ji,jj) = ( zt / z1_T0 + zn / zd ) * ztm
1168            !
1169      END_2D
1170      !
1171      IF( ln_timing )   CALL timing_stop('eos_pt_from_ct')
1172      !
1173   END FUNCTION eos_pt_from_ct
1174
1175
1176   SUBROUTINE eos_fzp_2d( psal, ptf, pdep )
1177      !!
1178      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   )           ::   psal   ! salinity   [psu]
1179      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ), OPTIONAL ::   pdep   ! depth      [m]
1180      REAL(wp), DIMENSION(:,:)    , INTENT(out  )           ::   ptf    ! freezing temperature [Celsius]
1181      !!
1182      CALL eos_fzp_2d_t( psal, ptf, is_tile(ptf), pdep )
1183   END SUBROUTINE eos_fzp_2d
1184
1185
1186   SUBROUTINE  eos_fzp_2d_t( psal, ptf, kttf, pdep )
1187      !!----------------------------------------------------------------------
1188      !!                 ***  ROUTINE eos_fzp  ***
1189      !!
1190      !! ** Purpose :   Compute the freezing point temperature [Celsius]
1191      !!
1192      !! ** Method  :   UNESCO freezing point (ptf) in Celsius is given by
1193      !!       ptf(t,z) = (-.0575+1.710523e-3*sqrt(abs(s))-2.154996e-4*s)*s - 7.53e-4*z
1194      !!       checkvalue: tf=-2.588567 Celsius for s=40psu, z=500m
1195      !!
1196      !! Reference  :   UNESCO tech. papers in the marine science no. 28. 1978
1197      !!----------------------------------------------------------------------
1198      INTEGER                       , INTENT(in   )           ::   kttf
1199      REAL(wp), DIMENSION(jpi,jpj)  , INTENT(in   )           ::   psal   ! salinity   [psu]
1200      REAL(wp), DIMENSION(jpi,jpj)  , INTENT(in   ), OPTIONAL ::   pdep   ! depth      [m]
1201      REAL(wp), DIMENSION(A2D_T(kttf)), INTENT(out  )           ::   ptf    ! freezing temperature [Celsius]
1202      !
1203      INTEGER  ::   ji, jj          ! dummy loop indices
1204      REAL(wp) ::   zt, zs, z1_S0   ! local scalars
1205      !!----------------------------------------------------------------------
1206      !
1207      SELECT CASE ( neos )
1208      !
1209      CASE ( np_teos10, np_seos )      !==  CT,SA (TEOS-10 and S-EOS formulations) ==!
1210         !
1211         z1_S0 = 1._wp / 35.16504_wp
1212         DO_2D( 1, 1, 1, 1 )
1213            zs= SQRT( ABS( psal(ji,jj) ) * z1_S0 )           ! square root salinity
1214            ptf(ji,jj) = ((((1.46873e-03_wp*zs-9.64972e-03_wp)*zs+2.28348e-02_wp)*zs &
1215               &          - 3.12775e-02_wp)*zs+2.07679e-02_wp)*zs-5.87701e-02_wp
1216         END_2D
1217         ptf(:,:) = ptf(:,:) * psal(:,:)
1218         !
1219         IF( PRESENT( pdep ) )   ptf(:,:) = ptf(:,:) - 7.53e-4 * pdep(:,:)
1220         !
1221      CASE ( np_eos80 )                !==  PT,SP (UNESCO formulation)  ==!
1222         !
1223         ptf(:,:) = ( - 0.0575_wp + 1.710523e-3_wp * SQRT( psal(:,:) )   &
1224            &                     - 2.154996e-4_wp *       psal(:,:)   ) * psal(:,:)
1225            !
1226         IF( PRESENT( pdep ) )   ptf(:,:) = ptf(:,:) - 7.53e-4 * pdep(:,:)
1227         !
1228      CASE DEFAULT
1229         WRITE(ctmp1,*) '          bad flag value for neos = ', neos
1230         CALL ctl_stop( 'eos_fzp_2d:', ctmp1 )
1231         !
1232      END SELECT
1233      !
1234  END SUBROUTINE eos_fzp_2d_t
1235
1236
1237  SUBROUTINE eos_fzp_0d( psal, ptf, pdep )
1238      !!----------------------------------------------------------------------
1239      !!                 ***  ROUTINE eos_fzp  ***
1240      !!
1241      !! ** Purpose :   Compute the freezing point temperature [Celsius]
1242      !!
1243      !! ** Method  :   UNESCO freezing point (ptf) in Celsius is given by
1244      !!       ptf(t,z) = (-.0575+1.710523e-3*sqrt(abs(s))-2.154996e-4*s)*s - 7.53e-4*z
1245      !!       checkvalue: tf=-2.588567 Celsius for s=40psu, z=500m
1246      !!
1247      !! Reference  :   UNESCO tech. papers in the marine science no. 28. 1978
1248      !!----------------------------------------------------------------------
1249      REAL(wp), INTENT(in )           ::   psal         ! salinity   [psu]
1250      REAL(wp), INTENT(in ), OPTIONAL ::   pdep         ! depth      [m]
1251      REAL(wp), INTENT(out)           ::   ptf          ! freezing temperature [Celsius]
1252      !
1253      REAL(wp) :: zs   ! local scalars
1254      !!----------------------------------------------------------------------
1255      !
1256      SELECT CASE ( neos )
1257      !
1258      CASE ( np_teos10, np_seos )      !==  CT,SA (TEOS-10 and S-EOS formulations) ==!
1259         !
1260         zs  = SQRT( ABS( psal ) / 35.16504_wp )           ! square root salinity
1261         ptf = ((((1.46873e-03_wp*zs-9.64972e-03_wp)*zs+2.28348e-02_wp)*zs &
1262                  &          - 3.12775e-02_wp)*zs+2.07679e-02_wp)*zs-5.87701e-02_wp
1263         ptf = ptf * psal
1264         !
1265         IF( PRESENT( pdep ) )   ptf = ptf - 7.53e-4 * pdep
1266         !
1267      CASE ( np_eos80 )                !==  PT,SP (UNESCO formulation)  ==!
1268         !
1269         ptf = ( - 0.0575_wp + 1.710523e-3_wp * SQRT( psal )   &
1270            &                - 2.154996e-4_wp *       psal   ) * psal
1271            !
1272         IF( PRESENT( pdep ) )   ptf = ptf - 7.53e-4 * pdep
1273         !
1274      CASE DEFAULT
1275         WRITE(ctmp1,*) '          bad flag value for neos = ', neos
1276         CALL ctl_stop( 'eos_fzp_0d:', ctmp1 )
1277         !
1278      END SELECT
1279      !
1280   END SUBROUTINE eos_fzp_0d
1281
1282
1283   SUBROUTINE eos_pen( pts, pab_pe, ppen, Kmm )
1284      !!----------------------------------------------------------------------
1285      !!                 ***  ROUTINE eos_pen  ***
1286      !!
1287      !! ** Purpose :   Calculates nonlinear anomalies of alpha_PE, beta_PE and PE at T-points
1288      !!
1289      !! ** Method  :   PE is defined analytically as the vertical
1290      !!                   primitive of EOS times -g integrated between 0 and z>0.
1291      !!                pen is the nonlinear bsq-PE anomaly: pen = ( PE - rho0 gz ) / rho0 gz - rd
1292      !!                                                      = 1/z * /int_0^z rd dz - rd
1293      !!                                where rd is the density anomaly (see eos_rhd function)
1294      !!                ab_pe are partial derivatives of PE anomaly with respect to T and S:
1295      !!                    ab_pe(1) = - 1/(rho0 gz) * dPE/dT + drd/dT = - d(pen)/dT
1296      !!                    ab_pe(2) =   1/(rho0 gz) * dPE/dS + drd/dS =   d(pen)/dS
1297      !!
1298      !! ** Action  : - pen         : PE anomaly given at T-points
1299      !!            : - pab_pe  : given at T-points
1300      !!                    pab_pe(:,:,:,jp_tem) is alpha_pe
1301      !!                    pab_pe(:,:,:,jp_sal) is beta_pe
1302      !!----------------------------------------------------------------------
1303      INTEGER                              , INTENT(in   ) ::   Kmm   ! time level index
1304      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in   ) ::   pts     ! pot. temperature & salinity
1305      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(  out) ::   pab_pe  ! alpha_pe and beta_pe
1306      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(  out) ::   ppen     ! potential energy anomaly
1307      !
1308      INTEGER  ::   ji, jj, jk                ! dummy loop indices
1309      REAL(wp) ::   zt , zh , zs , ztm        ! local scalars
1310      REAL(wp) ::   zn , zn0, zn1, zn2        !   -      -
1311      !!----------------------------------------------------------------------
1312      !
1313      IF( ln_timing )   CALL timing_start('eos_pen')
1314      !
1315      SELECT CASE ( neos )
1316      !
1317      CASE( np_teos10, np_eos80 )                !==  polynomial TEOS-10 / EOS-80 ==!
1318         !
1319         DO_3D( 1, 1, 1, 1, 1, jpkm1 )
1320            !
1321            zh  = gdept(ji,jj,jk,Kmm) * r1_Z0                                ! depth
1322            zt  = pts (ji,jj,jk,jp_tem) * r1_T0                           ! temperature
1323            zs  = SQRT( ABS( pts(ji,jj,jk,jp_sal) + rdeltaS ) * r1_S0 )   ! square root salinity
1324            ztm = tmask(ji,jj,jk)                                         ! tmask
1325            !
1326            ! potential energy non-linear anomaly
1327            zn2 = (PEN012)*zt   &
1328               &   + PEN102*zs+PEN002
1329               !
1330            zn1 = ((PEN021)*zt   &
1331               &   + PEN111*zs+PEN011)*zt   &
1332               &   + (PEN201*zs+PEN101)*zs+PEN001
1333               !
1334            zn0 = ((((PEN040)*zt   &
1335               &   + PEN130*zs+PEN030)*zt   &
1336               &   + (PEN220*zs+PEN120)*zs+PEN020)*zt   &
1337               &   + ((PEN310*zs+PEN210)*zs+PEN110)*zs+PEN010)*zt   &
1338               &   + (((PEN400*zs+PEN300)*zs+PEN200)*zs+PEN100)*zs+PEN000
1339               !
1340            zn  = ( zn2 * zh + zn1 ) * zh + zn0
1341            !
1342            ppen(ji,jj,jk)  = zn * zh * r1_rho0 * ztm
1343            !
1344            ! alphaPE non-linear anomaly
1345            zn2 = APE002
1346            !
1347            zn1 = (APE011)*zt   &
1348               &   + APE101*zs+APE001
1349               !
1350            zn0 = (((APE030)*zt   &
1351               &   + APE120*zs+APE020)*zt   &
1352               &   + (APE210*zs+APE110)*zs+APE010)*zt   &
1353               &   + ((APE300*zs+APE200)*zs+APE100)*zs+APE000
1354               !
1355            zn  = ( zn2 * zh + zn1 ) * zh + zn0
1356            !
1357            pab_pe(ji,jj,jk,jp_tem) = zn * zh * r1_rho0 * ztm
1358            !
1359            ! betaPE non-linear anomaly
1360            zn2 = BPE002
1361            !
1362            zn1 = (BPE011)*zt   &
1363               &   + BPE101*zs+BPE001
1364               !
1365            zn0 = (((BPE030)*zt   &
1366               &   + BPE120*zs+BPE020)*zt   &
1367               &   + (BPE210*zs+BPE110)*zs+BPE010)*zt   &
1368               &   + ((BPE300*zs+BPE200)*zs+BPE100)*zs+BPE000
1369               !
1370            zn  = ( zn2 * zh + zn1 ) * zh + zn0
1371            !
1372            pab_pe(ji,jj,jk,jp_sal) = zn / zs * zh * r1_rho0 * ztm
1373            !
1374         END_3D
1375         !
1376      CASE( np_seos )                !==  Vallis (2006) simplified EOS  ==!
1377         !
1378         DO_3D( 1, 1, 1, 1, 1, jpkm1 )
1379            zt  = pts(ji,jj,jk,jp_tem) - 10._wp  ! temperature anomaly (t-T0)
1380            zs = pts (ji,jj,jk,jp_sal) - 35._wp  ! abs. salinity anomaly (s-S0)
1381            zh  = gdept(ji,jj,jk,Kmm)              ! depth in meters  at t-point
1382            ztm = tmask(ji,jj,jk)                ! tmask
1383            zn  = 0.5_wp * zh * r1_rho0 * ztm
1384            !                                    ! Potential Energy
1385            ppen(ji,jj,jk) = ( rn_a0 * rn_mu1 * zt + rn_b0 * rn_mu2 * zs ) * zn
1386            !                                    ! alphaPE
1387            pab_pe(ji,jj,jk,jp_tem) = - rn_a0 * rn_mu1 * zn
1388            pab_pe(ji,jj,jk,jp_sal) =   rn_b0 * rn_mu2 * zn
1389            !
1390         END_3D
1391         !
1392      CASE( np_eeos )                !==  exponential SEAMOUNT EOS  ==!
1393         !
1394         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
1395            ppen(ji,jk,jk) = 0.
1396            !                                 
1397            pab_pe(ji,jj,jk,jp_tem) = 0.
1398            pab_pe(ji,jj,jk,jp_sal) = 0.
1399         END_3D
1400         !
1401      CASE DEFAULT
1402         WRITE(ctmp1,*) '          bad flag value for neos = ', neos
1403         CALL ctl_stop( 'eos_pen:', ctmp1 )
1404         !
1405      END SELECT
1406      !
1407      IF( ln_timing )   CALL timing_stop('eos_pen')
1408      !
1409   END SUBROUTINE eos_pen
1410
1411
1412   SUBROUTINE eos_init
1413      !!----------------------------------------------------------------------
1414      !!                 ***  ROUTINE eos_init  ***
1415      !!
1416      !! ** Purpose :   initializations for the equation of state
1417      !!
1418      !! ** Method  :   Read the namelist nameos and control the parameters
1419      !!----------------------------------------------------------------------
1420      INTEGER  ::   ios   ! local integer
1421      INTEGER  ::   ioptio   ! local integer
1422      !!
1423      NAMELIST/nameos/ ln_TEOS10, ln_EOS80, ln_SEOS, ln_EEOS, rn_a0, rn_b0, rn_lambda1, rn_mu1,   &
1424         &                                             rn_lambda2, rn_mu2, rn_nu
1425      !!----------------------------------------------------------------------
1426      !
1427      READ  ( numnam_ref, nameos, IOSTAT = ios, ERR = 901 )
1428901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nameos in reference namelist' )
1429      !
1430      READ  ( numnam_cfg, nameos, IOSTAT = ios, ERR = 902 )
1431902   IF( ios >  0 )   CALL ctl_nam ( ios , 'nameos in configuration namelist' )
1432      IF(lwm) WRITE( numond, nameos )
1433      !
1434      rho0        = 1000._wp                 !: volumic mass of reference     [kg/m3]
1435      rcp         = 3991.86795711963_wp      !: heat capacity     [J/K]
1436      !
1437      IF(lwp) THEN                ! Control print
1438         WRITE(numout,*)
1439         WRITE(numout,*) 'eos_init : equation of state'
1440         WRITE(numout,*) '~~~~~~~~'
1441         WRITE(numout,*) '   Namelist nameos : Chosen the Equation Of Seawater (EOS)'
1442         WRITE(numout,*) '      TEOS-10 : rho=F(Conservative Temperature, Absolute  Salinity, depth)   ln_TEOS10 = ', ln_TEOS10
1443         WRITE(numout,*) '      EOS-80  : rho=F(Potential    Temperature, Practical Salinity, depth)   ln_EOS80  = ', ln_EOS80
1444         WRITE(numout,*) '      S-EOS   : rho=F(Conservative Temperature, Absolute  Salinity, depth)   ln_SEOS   = ', ln_SEOS
1445         WRITE(numout,*) '      E-EOS   : rho=F(Conservative Temperature,                     depth)   ln_EEOS   = ', ln_EEOS
1446      ENDIF
1447
1448      ! Check options for equation of state & set neos based on logical flags
1449      ioptio = 0
1450      IF( ln_TEOS10 ) THEN   ;   ioptio = ioptio+1   ;   neos = np_teos10   ;   ENDIF
1451      IF( ln_EOS80  ) THEN   ;   ioptio = ioptio+1   ;   neos = np_eos80    ;   ENDIF
1452      IF( ln_SEOS   ) THEN   ;   ioptio = ioptio+1   ;   neos = np_seos     ;   ENDIF
1453      IF( ln_EEOS   ) THEN   ;   ioptio = ioptio+1   ;   neos = np_eeos     ;   ENDIF
1454      IF( ioptio /= 1 )   CALL ctl_stop("Exactly one equation of state option must be selected")
1455      !
1456      SELECT CASE( neos )         ! check option
1457      !
1458      CASE( np_teos10 )                       !==  polynomial TEOS-10  ==!
1459         IF(lwp) WRITE(numout,*)
1460         IF(lwp) WRITE(numout,*) '   ==>>>   use of TEOS-10 equation of state (cons. temp. and abs. salinity)'
1461         !
1462         l_useCT = .TRUE.                          ! model temperature is Conservative temperature
1463         !
1464         rdeltaS = 32._wp
1465         r1_S0  = 0.875_wp/35.16504_wp
1466         r1_T0  = 1._wp/40._wp
1467         r1_Z0  = 1.e-4_wp
1468         !
1469         EOS000 = 8.0189615746e+02_wp
1470         EOS100 = 8.6672408165e+02_wp
1471         EOS200 = -1.7864682637e+03_wp
1472         EOS300 = 2.0375295546e+03_wp
1473         EOS400 = -1.2849161071e+03_wp
1474         EOS500 = 4.3227585684e+02_wp
1475         EOS600 = -6.0579916612e+01_wp
1476         EOS010 = 2.6010145068e+01_wp
1477         EOS110 = -6.5281885265e+01_wp
1478         EOS210 = 8.1770425108e+01_wp
1479         EOS310 = -5.6888046321e+01_wp
1480         EOS410 = 1.7681814114e+01_wp
1481         EOS510 = -1.9193502195_wp
1482         EOS020 = -3.7074170417e+01_wp
1483         EOS120 = 6.1548258127e+01_wp
1484         EOS220 = -6.0362551501e+01_wp
1485         EOS320 = 2.9130021253e+01_wp
1486         EOS420 = -5.4723692739_wp
1487         EOS030 = 2.1661789529e+01_wp
1488         EOS130 = -3.3449108469e+01_wp
1489         EOS230 = 1.9717078466e+01_wp
1490         EOS330 = -3.1742946532_wp
1491         EOS040 = -8.3627885467_wp
1492         EOS140 = 1.1311538584e+01_wp
1493         EOS240 = -5.3563304045_wp
1494         EOS050 = 5.4048723791e-01_wp
1495         EOS150 = 4.8169980163e-01_wp
1496         EOS060 = -1.9083568888e-01_wp
1497         EOS001 = 1.9681925209e+01_wp
1498         EOS101 = -4.2549998214e+01_wp
1499         EOS201 = 5.0774768218e+01_wp
1500         EOS301 = -3.0938076334e+01_wp
1501         EOS401 = 6.6051753097_wp
1502         EOS011 = -1.3336301113e+01_wp
1503         EOS111 = -4.4870114575_wp
1504         EOS211 = 5.0042598061_wp
1505         EOS311 = -6.5399043664e-01_wp
1506         EOS021 = 6.7080479603_wp
1507         EOS121 = 3.5063081279_wp
1508         EOS221 = -1.8795372996_wp
1509         EOS031 = -2.4649669534_wp
1510         EOS131 = -5.5077101279e-01_wp
1511         EOS041 = 5.5927935970e-01_wp
1512         EOS002 = 2.0660924175_wp
1513         EOS102 = -4.9527603989_wp
1514         EOS202 = 2.5019633244_wp
1515         EOS012 = 2.0564311499_wp
1516         EOS112 = -2.1311365518e-01_wp
1517         EOS022 = -1.2419983026_wp
1518         EOS003 = -2.3342758797e-02_wp
1519         EOS103 = -1.8507636718e-02_wp
1520         EOS013 = 3.7969820455e-01_wp
1521         !
1522         ALP000 = -6.5025362670e-01_wp
1523         ALP100 = 1.6320471316_wp
1524         ALP200 = -2.0442606277_wp
1525         ALP300 = 1.4222011580_wp
1526         ALP400 = -4.4204535284e-01_wp
1527         ALP500 = 4.7983755487e-02_wp
1528         ALP010 = 1.8537085209_wp
1529         ALP110 = -3.0774129064_wp
1530         ALP210 = 3.0181275751_wp
1531         ALP310 = -1.4565010626_wp
1532         ALP410 = 2.7361846370e-01_wp
1533         ALP020 = -1.6246342147_wp
1534         ALP120 = 2.5086831352_wp
1535         ALP220 = -1.4787808849_wp
1536         ALP320 = 2.3807209899e-01_wp
1537         ALP030 = 8.3627885467e-01_wp
1538         ALP130 = -1.1311538584_wp
1539         ALP230 = 5.3563304045e-01_wp
1540         ALP040 = -6.7560904739e-02_wp
1541         ALP140 = -6.0212475204e-02_wp
1542         ALP050 = 2.8625353333e-02_wp
1543         ALP001 = 3.3340752782e-01_wp
1544         ALP101 = 1.1217528644e-01_wp
1545         ALP201 = -1.2510649515e-01_wp
1546         ALP301 = 1.6349760916e-02_wp
1547         ALP011 = -3.3540239802e-01_wp
1548         ALP111 = -1.7531540640e-01_wp
1549         ALP211 = 9.3976864981e-02_wp
1550         ALP021 = 1.8487252150e-01_wp
1551         ALP121 = 4.1307825959e-02_wp
1552         ALP031 = -5.5927935970e-02_wp
1553         ALP002 = -5.1410778748e-02_wp
1554         ALP102 = 5.3278413794e-03_wp
1555         ALP012 = 6.2099915132e-02_wp
1556         ALP003 = -9.4924551138e-03_wp
1557         !
1558         BET000 = 1.0783203594e+01_wp
1559         BET100 = -4.4452095908e+01_wp
1560         BET200 = 7.6048755820e+01_wp
1561         BET300 = -6.3944280668e+01_wp
1562         BET400 = 2.6890441098e+01_wp
1563         BET500 = -4.5221697773_wp
1564         BET010 = -8.1219372432e-01_wp
1565         BET110 = 2.0346663041_wp
1566         BET210 = -2.1232895170_wp
1567         BET310 = 8.7994140485e-01_wp
1568         BET410 = -1.1939638360e-01_wp
1569         BET020 = 7.6574242289e-01_wp
1570         BET120 = -1.5019813020_wp
1571         BET220 = 1.0872489522_wp
1572         BET320 = -2.7233429080e-01_wp
1573         BET030 = -4.1615152308e-01_wp
1574         BET130 = 4.9061350869e-01_wp
1575         BET230 = -1.1847737788e-01_wp
1576         BET040 = 1.4073062708e-01_wp
1577         BET140 = -1.3327978879e-01_wp
1578         BET050 = 5.9929880134e-03_wp
1579         BET001 = -5.2937873009e-01_wp
1580         BET101 = 1.2634116779_wp
1581         BET201 = -1.1547328025_wp
1582         BET301 = 3.2870876279e-01_wp
1583         BET011 = -5.5824407214e-02_wp
1584         BET111 = 1.2451933313e-01_wp
1585         BET211 = -2.4409539932e-02_wp
1586         BET021 = 4.3623149752e-02_wp
1587         BET121 = -4.6767901790e-02_wp
1588         BET031 = -6.8523260060e-03_wp
1589         BET002 = -6.1618945251e-02_wp
1590         BET102 = 6.2255521644e-02_wp
1591         BET012 = -2.6514181169e-03_wp
1592         BET003 = -2.3025968587e-04_wp
1593         !
1594         PEN000 = -9.8409626043_wp
1595         PEN100 = 2.1274999107e+01_wp
1596         PEN200 = -2.5387384109e+01_wp
1597         PEN300 = 1.5469038167e+01_wp
1598         PEN400 = -3.3025876549_wp
1599         PEN010 = 6.6681505563_wp
1600         PEN110 = 2.2435057288_wp
1601         PEN210 = -2.5021299030_wp
1602         PEN310 = 3.2699521832e-01_wp
1603         PEN020 = -3.3540239802_wp
1604         PEN120 = -1.7531540640_wp
1605         PEN220 = 9.3976864981e-01_wp
1606         PEN030 = 1.2324834767_wp
1607         PEN130 = 2.7538550639e-01_wp
1608         PEN040 = -2.7963967985e-01_wp
1609         PEN001 = -1.3773949450_wp
1610         PEN101 = 3.3018402659_wp
1611         PEN201 = -1.6679755496_wp
1612         PEN011 = -1.3709540999_wp
1613         PEN111 = 1.4207577012e-01_wp
1614         PEN021 = 8.2799886843e-01_wp
1615         PEN002 = 1.7507069098e-02_wp
1616         PEN102 = 1.3880727538e-02_wp
1617         PEN012 = -2.8477365341e-01_wp
1618         !
1619         APE000 = -1.6670376391e-01_wp
1620         APE100 = -5.6087643219e-02_wp
1621         APE200 = 6.2553247576e-02_wp
1622         APE300 = -8.1748804580e-03_wp
1623         APE010 = 1.6770119901e-01_wp
1624         APE110 = 8.7657703198e-02_wp
1625         APE210 = -4.6988432490e-02_wp
1626         APE020 = -9.2436260751e-02_wp
1627         APE120 = -2.0653912979e-02_wp
1628         APE030 = 2.7963967985e-02_wp
1629         APE001 = 3.4273852498e-02_wp
1630         APE101 = -3.5518942529e-03_wp
1631         APE011 = -4.1399943421e-02_wp
1632         APE002 = 7.1193413354e-03_wp
1633         !
1634         BPE000 = 2.6468936504e-01_wp
1635         BPE100 = -6.3170583896e-01_wp
1636         BPE200 = 5.7736640125e-01_wp
1637         BPE300 = -1.6435438140e-01_wp
1638         BPE010 = 2.7912203607e-02_wp
1639         BPE110 = -6.2259666565e-02_wp
1640         BPE210 = 1.2204769966e-02_wp
1641         BPE020 = -2.1811574876e-02_wp
1642         BPE120 = 2.3383950895e-02_wp
1643         BPE030 = 3.4261630030e-03_wp
1644         BPE001 = 4.1079296834e-02_wp
1645         BPE101 = -4.1503681096e-02_wp
1646         BPE011 = 1.7676120780e-03_wp
1647         BPE002 = 1.7269476440e-04_wp
1648         !
1649      CASE( np_eos80 )                        !==  polynomial EOS-80 formulation  ==!
1650         !
1651         IF(lwp) WRITE(numout,*)
1652         IF(lwp) WRITE(numout,*) '   ==>>>   use of EOS-80 equation of state (pot. temp. and pract. salinity)'
1653         !
1654         l_useCT = .FALSE.                         ! model temperature is Potential temperature
1655         rdeltaS = 20._wp
1656         r1_S0  = 1._wp/40._wp
1657         r1_T0  = 1._wp/40._wp
1658         r1_Z0  = 1.e-4_wp
1659         !
1660         EOS000 = 9.5356891948e+02_wp
1661         EOS100 = 1.7136499189e+02_wp
1662         EOS200 = -3.7501039454e+02_wp
1663         EOS300 = 5.1856810420e+02_wp
1664         EOS400 = -3.7264470465e+02_wp
1665         EOS500 = 1.4302533998e+02_wp
1666         EOS600 = -2.2856621162e+01_wp
1667         EOS010 = 1.0087518651e+01_wp
1668         EOS110 = -1.3647741861e+01_wp
1669         EOS210 = 8.8478359933_wp
1670         EOS310 = -7.2329388377_wp
1671         EOS410 = 1.4774410611_wp
1672         EOS510 = 2.0036720553e-01_wp
1673         EOS020 = -2.5579830599e+01_wp
1674         EOS120 = 2.4043512327e+01_wp
1675         EOS220 = -1.6807503990e+01_wp
1676         EOS320 = 8.3811577084_wp
1677         EOS420 = -1.9771060192_wp
1678         EOS030 = 1.6846451198e+01_wp
1679         EOS130 = -2.1482926901e+01_wp
1680         EOS230 = 1.0108954054e+01_wp
1681         EOS330 = -6.2675951440e-01_wp
1682         EOS040 = -8.0812310102_wp
1683         EOS140 = 1.0102374985e+01_wp
1684         EOS240 = -4.8340368631_wp
1685         EOS050 = 1.2079167803_wp
1686         EOS150 = 1.1515380987e-01_wp
1687         EOS060 = -2.4520288837e-01_wp
1688         EOS001 = 1.0748601068e+01_wp
1689         EOS101 = -1.7817043500e+01_wp
1690         EOS201 = 2.2181366768e+01_wp
1691         EOS301 = -1.6750916338e+01_wp
1692         EOS401 = 4.1202230403_wp
1693         EOS011 = -1.5852644587e+01_wp
1694         EOS111 = -7.6639383522e-01_wp
1695         EOS211 = 4.1144627302_wp
1696         EOS311 = -6.6955877448e-01_wp
1697         EOS021 = 9.9994861860_wp
1698         EOS121 = -1.9467067787e-01_wp
1699         EOS221 = -1.2177554330_wp
1700         EOS031 = -3.4866102017_wp
1701         EOS131 = 2.2229155620e-01_wp
1702         EOS041 = 5.9503008642e-01_wp
1703         EOS002 = 1.0375676547_wp
1704         EOS102 = -3.4249470629_wp
1705         EOS202 = 2.0542026429_wp
1706         EOS012 = 2.1836324814_wp
1707         EOS112 = -3.4453674320e-01_wp
1708         EOS022 = -1.2548163097_wp
1709         EOS003 = 1.8729078427e-02_wp
1710         EOS103 = -5.7238495240e-02_wp
1711         EOS013 = 3.8306136687e-01_wp
1712         !
1713         ALP000 = -2.5218796628e-01_wp
1714         ALP100 = 3.4119354654e-01_wp
1715         ALP200 = -2.2119589983e-01_wp
1716         ALP300 = 1.8082347094e-01_wp
1717         ALP400 = -3.6936026529e-02_wp
1718         ALP500 = -5.0091801383e-03_wp
1719         ALP010 = 1.2789915300_wp
1720         ALP110 = -1.2021756164_wp
1721         ALP210 = 8.4037519952e-01_wp
1722         ALP310 = -4.1905788542e-01_wp
1723         ALP410 = 9.8855300959e-02_wp
1724         ALP020 = -1.2634838399_wp
1725         ALP120 = 1.6112195176_wp
1726         ALP220 = -7.5817155402e-01_wp
1727         ALP320 = 4.7006963580e-02_wp
1728         ALP030 = 8.0812310102e-01_wp
1729         ALP130 = -1.0102374985_wp
1730         ALP230 = 4.8340368631e-01_wp
1731         ALP040 = -1.5098959754e-01_wp
1732         ALP140 = -1.4394226233e-02_wp
1733         ALP050 = 3.6780433255e-02_wp
1734         ALP001 = 3.9631611467e-01_wp
1735         ALP101 = 1.9159845880e-02_wp
1736         ALP201 = -1.0286156825e-01_wp
1737         ALP301 = 1.6738969362e-02_wp
1738         ALP011 = -4.9997430930e-01_wp
1739         ALP111 = 9.7335338937e-03_wp
1740         ALP211 = 6.0887771651e-02_wp
1741         ALP021 = 2.6149576513e-01_wp
1742         ALP121 = -1.6671866715e-02_wp
1743         ALP031 = -5.9503008642e-02_wp
1744         ALP002 = -5.4590812035e-02_wp
1745         ALP102 = 8.6134185799e-03_wp
1746         ALP012 = 6.2740815484e-02_wp
1747         ALP003 = -9.5765341718e-03_wp
1748         !
1749         BET000 = 2.1420623987_wp
1750         BET100 = -9.3752598635_wp
1751         BET200 = 1.9446303907e+01_wp
1752         BET300 = -1.8632235232e+01_wp
1753         BET400 = 8.9390837485_wp
1754         BET500 = -1.7142465871_wp
1755         BET010 = -1.7059677327e-01_wp
1756         BET110 = 2.2119589983e-01_wp
1757         BET210 = -2.7123520642e-01_wp
1758         BET310 = 7.3872053057e-02_wp
1759         BET410 = 1.2522950346e-02_wp
1760         BET020 = 3.0054390409e-01_wp
1761         BET120 = -4.2018759976e-01_wp
1762         BET220 = 3.1429341406e-01_wp
1763         BET320 = -9.8855300959e-02_wp
1764         BET030 = -2.6853658626e-01_wp
1765         BET130 = 2.5272385134e-01_wp
1766         BET230 = -2.3503481790e-02_wp
1767         BET040 = 1.2627968731e-01_wp
1768         BET140 = -1.2085092158e-01_wp
1769         BET050 = 1.4394226233e-03_wp
1770         BET001 = -2.2271304375e-01_wp
1771         BET101 = 5.5453416919e-01_wp
1772         BET201 = -6.2815936268e-01_wp
1773         BET301 = 2.0601115202e-01_wp
1774         BET011 = -9.5799229402e-03_wp
1775         BET111 = 1.0286156825e-01_wp
1776         BET211 = -2.5108454043e-02_wp
1777         BET021 = -2.4333834734e-03_wp
1778         BET121 = -3.0443885826e-02_wp
1779         BET031 = 2.7786444526e-03_wp
1780         BET002 = -4.2811838287e-02_wp
1781         BET102 = 5.1355066072e-02_wp
1782         BET012 = -4.3067092900e-03_wp
1783         BET003 = -7.1548119050e-04_wp
1784         !
1785         PEN000 = -5.3743005340_wp
1786         PEN100 = 8.9085217499_wp
1787         PEN200 = -1.1090683384e+01_wp
1788         PEN300 = 8.3754581690_wp
1789         PEN400 = -2.0601115202_wp
1790         PEN010 = 7.9263222935_wp
1791         PEN110 = 3.8319691761e-01_wp
1792         PEN210 = -2.0572313651_wp
1793         PEN310 = 3.3477938724e-01_wp
1794         PEN020 = -4.9997430930_wp
1795         PEN120 = 9.7335338937e-02_wp
1796         PEN220 = 6.0887771651e-01_wp
1797         PEN030 = 1.7433051009_wp
1798         PEN130 = -1.1114577810e-01_wp
1799         PEN040 = -2.9751504321e-01_wp
1800         PEN001 = -6.9171176978e-01_wp
1801         PEN101 = 2.2832980419_wp
1802         PEN201 = -1.3694684286_wp
1803         PEN011 = -1.4557549876_wp
1804         PEN111 = 2.2969116213e-01_wp
1805         PEN021 = 8.3654420645e-01_wp
1806         PEN002 = -1.4046808820e-02_wp
1807         PEN102 = 4.2928871430e-02_wp
1808         PEN012 = -2.8729602515e-01_wp
1809         !
1810         APE000 = -1.9815805734e-01_wp
1811         APE100 = -9.5799229402e-03_wp
1812         APE200 = 5.1430784127e-02_wp
1813         APE300 = -8.3694846809e-03_wp
1814         APE010 = 2.4998715465e-01_wp
1815         APE110 = -4.8667669469e-03_wp
1816         APE210 = -3.0443885826e-02_wp
1817         APE020 = -1.3074788257e-01_wp
1818         APE120 = 8.3359333577e-03_wp
1819         APE030 = 2.9751504321e-02_wp
1820         APE001 = 3.6393874690e-02_wp
1821         APE101 = -5.7422790533e-03_wp
1822         APE011 = -4.1827210323e-02_wp
1823         APE002 = 7.1824006288e-03_wp
1824         !
1825         BPE000 = 1.1135652187e-01_wp
1826         BPE100 = -2.7726708459e-01_wp
1827         BPE200 = 3.1407968134e-01_wp
1828         BPE300 = -1.0300557601e-01_wp
1829         BPE010 = 4.7899614701e-03_wp
1830         BPE110 = -5.1430784127e-02_wp
1831         BPE210 = 1.2554227021e-02_wp
1832         BPE020 = 1.2166917367e-03_wp
1833         BPE120 = 1.5221942913e-02_wp
1834         BPE030 = -1.3893222263e-03_wp
1835         BPE001 = 2.8541225524e-02_wp
1836         BPE101 = -3.4236710714e-02_wp
1837         BPE011 = 2.8711395266e-03_wp
1838         BPE002 = 5.3661089288e-04_wp
1839         !
1840      CASE( np_seos )                        !==  Simplified EOS     ==!
1841
1842         r1_S0  = 0.875_wp/35.16504_wp   ! Used to convert CT in potential temperature when using bulk formulae (eos_pt_from_ct)
1843
1844         IF(lwp) THEN
1845            WRITE(numout,*)
1846            WRITE(numout,*) '   ==>>>   use of simplified eos:    '
1847            WRITE(numout,*) '              rhd(dT=T-10,dS=S-35,Z) = [-a0*(1+lambda1/2*dT+mu1*Z)*dT '
1848            WRITE(numout,*) '                                       + b0*(1+lambda2/2*dT+mu2*Z)*dS - nu*dT*dS] / rho0'
1849            WRITE(numout,*) '              with the following coefficients :'
1850            WRITE(numout,*) '                 thermal exp. coef.    rn_a0      = ', rn_a0
1851            WRITE(numout,*) '                 saline  cont. coef.   rn_b0      = ', rn_b0
1852            WRITE(numout,*) '                 cabbeling coef.       rn_lambda1 = ', rn_lambda1
1853            WRITE(numout,*) '                 cabbeling coef.       rn_lambda2 = ', rn_lambda2
1854            WRITE(numout,*) '                 thermobar. coef.      rn_mu1     = ', rn_mu1
1855            WRITE(numout,*) '                 thermobar. coef.      rn_mu2     = ', rn_mu2
1856            WRITE(numout,*) '                 2nd cabbel. coef.     rn_nu      = ', rn_nu
1857            WRITE(numout,*) '              Caution: rn_beta0=0 incompatible with ddm parameterization '
1858         ENDIF
1859         l_useCT = .TRUE.          ! Use conservative temperature
1860         !
1861      CASE( np_eeos )                        !==  Exponential SEAMOUNT EOS     ==!
1862
1863         r1_S0  = 0.875_wp/35.16504_wp   ! Used to convert CT in potential temperature when using bulk formulae (eos_pt_from_ct)
1864         
1865         IF(lwp) THEN
1866            WRITE(numout,*)
1867            WRITE(numout,*) '          use of exponential SEAMOUNT eos:    rhd(dT=T-(exp_profile), Z) = '
1868            WRITE(numout,*) '             [ -a0*dT   ]/rho0'
1869            WRITE(numout,*)
1870            WRITE(numout,*) '             thermal exp. coef.    rn_a0      = ', rn_a0
1871         ENDIF
1872         l_useCT = .TRUE.          ! Use conservative temperature
1873         !
1874      CASE DEFAULT                     !==  ERROR in neos  ==!
1875         WRITE(ctmp1,*) '          bad flag value for neos = ', neos, '. You should never see this error'
1876         CALL ctl_stop( ctmp1 )
1877         !
1878      END SELECT
1879      !
1880      rho0_rcp    = rho0 * rcp
1881      r1_rho0     = 1._wp / rho0
1882      r1_rcp      = 1._wp / rcp
1883      r1_rho0_rcp = 1._wp / rho0_rcp
1884      !
1885      IF(lwp) THEN
1886         IF( l_useCT )   THEN
1887            WRITE(numout,*)
1888            WRITE(numout,*) '   ==>>>   model uses Conservative Temperature'
1889            WRITE(numout,*) '           Important: model must be initialized with CT and SA fields'
1890         ELSE
1891            WRITE(numout,*)
1892            WRITE(numout,*) '   ==>>>   model does not use Conservative Temperature'
1893         ENDIF
1894      ENDIF
1895      !
1896      IF(lwp) WRITE(numout,*)
1897      IF(lwp) WRITE(numout,*) '   Associated physical constant'
1898      IF(lwp) WRITE(numout,*) '      volumic mass of reference           rho0  = ', rho0   , ' kg/m^3'
1899      IF(lwp) WRITE(numout,*) '      1. / rho0                        r1_rho0  = ', r1_rho0, ' m^3/kg'
1900      IF(lwp) WRITE(numout,*) '      ocean specific heat                 rcp   = ', rcp    , ' J/Kelvin'
1901      IF(lwp) WRITE(numout,*) '      rho0 * rcp                       rho0_rcp = ', rho0_rcp
1902      IF(lwp) WRITE(numout,*) '      1. / ( rho0 * rcp )           r1_rho0_rcp = ', r1_rho0_rcp
1903      !
1904   END SUBROUTINE eos_init
1905
1906   !!======================================================================
1907END MODULE eosbn2
Note: See TracBrowser for help on using the repository browser.