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/2016/dev_r6409_SIMPLIF_2_usrdef_tools/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2016/dev_r6409_SIMPLIF_2_usrdef_tools/NEMOGCM/NEMO/OPA_SRC/TRA/eosbn2.F90 @ 6827

Last change on this file since 6827 was 6827, checked in by flavoni, 8 years ago

#1692 - branch SIMPLIF_2_usrdef: commit routines Fortran to create domain_cfg.nc file, to have domain (input) files for new simplification branch. To be moved in TOOLS directory

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