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/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/RK3_SRC/TRA – NEMO

source: branches/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/RK3_SRC/TRA/eosbn2.F90 @ 8568

Last change on this file since 8568 was 8568, checked in by gm, 7 years ago

#1911 (ENHANCE-09): PART I.2 - _NONE option + remove zts + see associated wiki page

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