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 branches/2015/dev_r5177_CNRS4_stopar/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2015/dev_r5177_CNRS4_stopar/NEMOGCM/NEMO/OPA_SRC/TRA/eosbn2.F90 @ 5296

Last change on this file since 5296 was 5296, checked in by pabouttier, 9 years ago

Commit stochastic parametrization module and perturbation of EOS

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