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 @ 14010

Last change on this file since 14010 was 14010, checked in by gsamson, 3 years ago

merge dev_r2052_ENHANCE-09_rbourdal_massfluxconvection into trunk

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