New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
eosbn2.F90 in NEMO/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE/TRA – NEMO

source: NEMO/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE/TRA/eosbn2.F90 @ 14286

Last change on this file since 14286 was 14219, checked in by mcastril, 3 years ago

Add Mixed Precision support by Oriol Tintó

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