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

source: NEMO/branches/2019/dev_r12072_MERGE_OPTION2_2019/tests/ISOMIP+/MY_SRC/eosbn2.F90 @ 12202

Last change on this file since 12202 was 12202, checked in by cetlod, 4 years ago

dev_merge_option2 : merge in dev_r11613_ENHANCE-04_namelists_as_internalfiles

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