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_r13541_TOP-01_rlod_Antarctic_ice_Sheet_Fe_Source/tests/ISOMIP+/MY_SRC – NEMO

source: NEMO/branches/2020/dev_r13541_TOP-01_rlod_Antarctic_ice_Sheet_Fe_Source/tests/ISOMIP+/MY_SRC/eosbn2.F90 @ 13886

Last change on this file since 13886 was 13886, checked in by rlod, 3 years ago

phasing with trunk at revision r13787

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