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

source: NEMO/branches/2020/dev_r12563_ASINTER-06_ABL_improvement/tests/ISOMIP+/MY_SRC/eosbn2.F90 @ 12587

Last change on this file since 12587 was 12489, checked in by davestorkey, 4 years ago

Preparation for new timestepping scheme #2390.
Main changes:

  1. Initial euler timestep now handled in stp and not in TRA/DYN routines.
  2. Renaming of all timestep parameters. In summary, the namelist parameter is now rn_Dt and the current timestep is rDt (and rDt_ice, rDt_trc etc).
  3. Renaming of a few miscellaneous parameters, eg. atfp -> rn_atfp (namelist parameter used everywhere) and rau0 -> rho0.

This version gives bit-comparable results to the previous version of the trunk.

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