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

source: branches/UKMO/dev_r5518_GO6_package_OMP/NEMOGCM/NEMO/OPA_SRC/TRA/eosbn2.F90 @ 9176

Last change on this file since 9176 was 9176, checked in by andmirek, 6 years ago

#2001: OMP directives

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