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/trunk/src/OCE/TRA – NEMO

source: NEMO/trunk/src/OCE/TRA/eosbn2.F90 @ 11989

Last change on this file since 11989 was 11989, checked in by cetlod, 4 years ago

cmip6 diag : finalisation of CMIP6 diags implementation

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