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

source: NEMO/trunk/tests/ISOMIP+/MY_SRC/eosbn2.F90

Last change on this file was 15145, checked in by smasson, 3 years ago

trunk: pass all sette tests in debug with nn_hls = 2

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