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 @ 4915

Last change on this file since 4915 was 4915, checked in by cetlod, 9 years ago

dev_r4615_CNRS1_10_TEOS10_Ediag : add missing lbc_lnk

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