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 branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/TRA/eosbn2.F90 @ 7037

Last change on this file since 7037 was 7037, checked in by mocavero, 8 years ago

ORCA2_LIM_PISCES hybrid version update

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