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

source: branches/2014/dev_CNRS1_10_TEOS10_Ediag/NEMOGCM/NEMO/OPA_SRC/TRA/eosbn2.F90 @ 4820

Last change on this file since 4820 was 4820, checked in by gm, 10 years ago

#1294 : TEOS-10 : correct wrong coefficients on beta following + style corrections

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