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 trunk/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: trunk/NEMOGCM/NEMO/OPA_SRC/TRA/eosbn2.F90 @ 6505

Last change on this file since 6505 was 6505, checked in by gm, 8 years ago

#1722 - trunk: TEOS-10 freezing point bug correction

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