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/r12581_ticket2418/tests/ISOMIP+/MY_SRC – NEMO

source: NEMO/branches/2020/r12581_ticket2418/tests/ISOMIP+/MY_SRC/eosbn2.F90 @ 12930

Last change on this file since 12930 was 12930, checked in by smasson, 5 years ago

r12581_ticket2418: update with trunk @12929, see #2418

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