New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
eosbn2.F90 in NEMO/branches/UKMO/dev_r12866_HPC-02_Daley_Tiling_trial_extra_halo/src/OCE/TRA – NEMO

source: NEMO/branches/UKMO/dev_r12866_HPC-02_Daley_Tiling_trial_extra_halo/src/OCE/TRA/eosbn2.F90 @ 13409

Last change on this file since 13409 was 13409, checked in by hadcv, 3 years ago

Remaining changes prior to trunk merge

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