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/UKMO/NEMO_4.0.2_GO8_package_ENHANCE-02_ISF_nemo/tests/ISOMIP+/MY_SRC – NEMO

source: NEMO/branches/UKMO/NEMO_4.0.2_GO8_package_ENHANCE-02_ISF_nemo/tests/ISOMIP+/MY_SRC/eosbn2.F90 @ 12721

Last change on this file since 12721 was 12721, checked in by mathiot, 4 years ago

NEMO_4.0.2_GO8_package_ENHANCE-02_ISF_nemo: add last year isf dev

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