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/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/TRA/eosbn2.F90 @ 4478

Last change on this file since 4478 was 4478, checked in by trackstand2, 10 years ago

Changes in eosbn2 to make precision of f.p. values explicit

  • Property svn:keywords set to Id
File size: 45.1 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
18   !!            3.3  ! 2010-05  (C. Ethe, G. Madec)  merge TRC-TRA
19   !!             -   ! 2010-10  (G. Nurser, G. Madec)  add eos_alpbet used in ldfslp
20   !!----------------------------------------------------------------------
21
22   !!----------------------------------------------------------------------
23   !!   eos            : generic interface of the equation of state
24   !!   eos_insitu     : Compute the in situ density
25   !!   eos_insitu_pot : Compute the insitu and surface referenced potential
26   !!                    volumic mass
27   !!   eos_insitu_2d  : Compute the in situ density for 2d fields
28   !!   eos_bn2        : Compute the Brunt-Vaisala frequency
29   !!   eos_alpbet     : calculates the in situ thermal and haline expansion coeff.
30   !!   tfreez         : Compute the surface freezing temperature
31   !!   eos_init       : set eos parameters (namelist)
32   !!----------------------------------------------------------------------
33   USE dom_oce         ! ocean space and time domain
34   USE phycst          ! physical constants
35   USE zdfddm          ! vertical physics: double diffusion
36   USE in_out_manager  ! I/O manager
37   USE lib_mpp         ! MPP library
38   USE prtctl          ! Print control
39
40   IMPLICIT NONE
41   PRIVATE
42
43   !                   !! * Interface
44   INTERFACE eos
45      MODULE PROCEDURE eos_insitu, eos_insitu_pot, eos_insitu_2d
46   END INTERFACE
47   INTERFACE bn2
48      MODULE PROCEDURE eos_bn2
49   END INTERFACE
50
51   PUBLIC   eos        ! called by step, istate, tranpc and zpsgrd modules
52   PUBLIC   eos_init   ! called by istate module
53   PUBLIC   bn2        ! called by step module
54   PUBLIC   eos_alpbet ! called by ldfslp module
55   PUBLIC   tfreez     ! called by sbcice_... modules
56
57   !                                          !!* Namelist (nameos) *
58   INTEGER , PUBLIC ::   nn_eos   = 0         !: = 0/1/2 type of eq. of state and Brunt-Vaisala frequ.
59   REAL(wp), PUBLIC ::   rn_alpha = 2.0e-4_wp !: thermal expension coeff. (linear equation of state)
60   REAL(wp), PUBLIC ::   rn_beta  = 7.7e-4_wp !: saline  expension coeff. (linear equation of state)
61
62   REAL(wp), PUBLIC ::   ralpbet              !: alpha / beta ratio
63
64   !! * Control permutation of array indices
65#  include "dom_oce_ftrans.h90"
66#  include "zdfddm_ftrans.h90"
67   
68   !! * Substitutions
69#  include "domzgr_substitute.h90"
70#  include "vectopt_loop_substitute.h90"
71   !!----------------------------------------------------------------------
72   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
73   !! $Id$
74   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
75   !!----------------------------------------------------------------------
76CONTAINS
77
78   SUBROUTINE eos_insitu( pts, prd )
79      !!----------------------------------------------------------------------
80      !!                   ***  ROUTINE eos_insitu  ***
81      !!
82      !! ** Purpose :   Compute the in situ density (ratio rho/rau0) from
83      !!       potential temperature and salinity using an equation of state
84      !!       defined through the namelist parameter nn_eos.
85      !!
86      !! ** Method  :   3 cases:
87      !!      nn_eos = 0 : Jackett and McDougall (1994) equation of state.
88      !!         the in situ density is computed directly as a function of
89      !!         potential temperature relative to the surface (the opa t
90      !!         variable), salt and pressure (assuming no pressure variation
91      !!         along geopotential surfaces, i.e. the pressure p in decibars
92      !!         is approximated by the depth in meters.
93      !!              prd(t,s,p) = ( rho(t,s,p) - rau0 ) / rau0
94      !!         with pressure                      p        decibars
95      !!              potential temperature         t        deg celsius
96      !!              salinity                      s        psu
97      !!              reference volumic mass        rau0     kg/m**3
98      !!              in situ volumic mass          rho      kg/m**3
99      !!              in situ density anomalie      prd      no units
100      !!         Check value: rho = 1060.93298 kg/m**3 for p=10000 dbar,
101      !!          t = 40 deg celcius, s=40 psu
102      !!      nn_eos = 1 : linear equation of state function of temperature only
103      !!              prd(t) = 0.0285 - rn_alpha * t
104      !!      nn_eos = 2 : linear equation of state function of temperature and
105      !!               salinity
106      !!              prd(t,s) = rn_beta * s - rn_alpha * tn - 1.
107      !!      Note that no boundary condition problem occurs in this routine
108      !!      as pts are defined over the whole domain.
109      !!
110      !! ** Action  :   compute prd , the in situ density (no units)
111      !!
112      !! References :   Jackett and McDougall, J. Atmos. Ocean. Tech., 1994
113      !!----------------------------------------------------------------------
114      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
115      USE wrk_nemo, ONLY:   zws => wrk_3d_1   ! 3D workspace
116      !!
117
118!FTRANS zws :I :I :z
119!FTRANS pts :I :I :z :I
120!FTRANS prd :I :I :z
121
122      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pts   ! 1 : potential temperature  [Celcius]
123      !                                                      ! 2 : salinity               [psu]
124      REAL(wp), DIMENSION(:,:,:)  , INTENT(  out) ::   prd   ! in situ density            [-]
125      !!
126      INTEGER  ::   ji, jj, jk           ! dummy loop indices
127      REAL(wp) ::   zt , zs , zh , zsr   ! local scalars
128      REAL(wp) ::   zr1, zr2, zr3, zr4   !   -      -
129      REAL(wp) ::   zrhop, ze, zbw, zb   !   -      -
130      REAL(wp) ::   zd , zc , zaw, za    !   -      -
131      REAL(wp) ::   zb1, za1, zkw, zk0   !   -      -
132      REAL(wp) ::   zrau0r               !   -      -
133      !!----------------------------------------------------------------------
134
135      IF( wrk_in_use(3, 1) ) THEN
136         CALL ctl_stop('eos_insitu: requested workspace array unavailable')   ;   RETURN
137      ENDIF
138
139      SELECT CASE( nn_eos )
140      !
141      CASE( 0 )                !==  Jackett and McDougall (1994) formulation  ==!
142         zrau0r = 1.e0 / rau0
143!CDIR NOVERRCHK
144#if defined key_z_first
145         DO jj = 1, jpj
146            DO ji = 1, jpi
147               DO jk = 1, mbkmax(ji,jj)
148                  zws(ji,jj,jk) = SQRT( ABS( pts(ji,jj,jk,jp_sal) ) )
149               END DO
150            END DO
151         END DO
152#else
153         zws(:,:,:) = SQRT( ABS( pts(:,:,:,jp_sal) ) )
154#endif
155         
156#if defined key_z_first
157         DO jj = 1, jpj
158            DO ji = 1, jpi
159               DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1
160#else
161         DO jk = 1, jpkm1
162            DO jj = 1, jpj
163               DO ji = 1, jpi
164#endif
165                  zt = pts   (ji,jj,jk,jp_tem)
166                  zs = pts   (ji,jj,jk,jp_sal)
167                  zh = fsdept(ji,jj,jk)        ! depth
168                  zsr= zws   (ji,jj,jk)        ! square root salinity
169                  !
170                  ! compute volumic mass pure water at atm pressure
171                  zr1= ( ( ( ( 6.536332e-9_wp  *zt - 1.120083e-6_wp )*zt + 1.001685e-4_wp )*zt   &
172                     &        -9.095290e-3_wp )*zt + 6.793952e-2_wp )*zt +  999.842594_wp
173                  ! seawater volumic mass atm pressure
174                  zr2= ( ( ( 5.3875e-9_wp*zt-8.2467e-7_wp ) *zt+7.6438e-5_wp ) *zt        &
175                     &                      -4.0899e-3_wp ) *zt+0.824493_wp
176                  zr3= ( -1.6546e-6_wp*zt+1.0227e-4_wp )    *zt-5.72466e-3_wp
177                  zr4= 4.8314e-4_wp
178                  !
179                  ! potential volumic mass (reference to the surface)
180                  zrhop= ( zr4*zs + zr3*zsr + zr2 ) *zs + zr1
181                  !
182                  ! add the compression terms
183                  ze = ( -3.508914e-8_wp*zt-1.248266e-8_wp ) *zt-2.595994e-6_wp
184                  zbw= (  1.296821e-6_wp*zt-5.782165e-9_wp ) *zt+1.045941e-4_wp
185                  zb = zbw + ze * zs
186                  !
187                  zd = -2.042967e-2_wp
188                  zc =   (-7.267926e-5_wp*zt+2.598241e-3_wp ) *zt+0.1571896_wp
189                  zaw= ( ( 5.939910e-6_wp*zt+2.512549e-3_wp ) *zt-0.1028859_wp ) *zt - 4.721788_wp
190                  za = ( zd*zsr + zc ) *zs + zaw
191                  !
192                  zb1=   (-0.1909078_wp*zt+7.390729_wp )        *zt-55.87545_wp
193                  za1= ( ( 2.326469e-3_wp*zt+1.553190_wp)       *zt-65.00517_wp ) *zt+1044.077_wp
194                  zkw= ( ( (-1.361629e-4_wp*zt-1.852732e-2_wp ) *zt-30.41638_wp ) *zt + 2098.925_wp ) *zt+190925.6_wp
195                  zk0= ( zb1*zsr + za1 )*zs + zkw
196                  !
197                  ! masked in situ density anomaly
198                  prd(ji,jj,jk) = (  zrhop / (  1.0_wp - zh / ( zk0 - zh * ( za - zh * zb ) )  )    &
199                     &             - rau0  ) * zrau0r * tmask(ji,jj,jk)
200               END DO
201            END DO
202         END DO
203         !
204      CASE( 1 )                !==  Linear formulation function of temperature only  ==!
205#if defined key_z_first
206         DO jj = 1, jpj
207            DO ji = 1, jpi
208               DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1
209                  prd(ji,jj,jk) = ( 0.0285_wp - rn_alpha * pts(ji,jj,jk,jp_tem) ) * tmask(ji,jj,jk)
210               END DO
211            END DO
212         END DO
213#else
214         DO jk = 1, jpkm1
215            prd(:,:,jk) = ( 0.0285_wp - rn_alpha * pts(:,:,jk,jp_tem) ) * tmask(:,:,jk)
216         END DO
217#endif
218         !
219      CASE( 2 )                !==  Linear formulation function of temperature and salinity  ==!
220#if defined key_z_first
221         DO jj = 1, jpj
222            DO ji = 1, jpi
223               DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1
224                  prd(ji,jj,jk) = ( rn_beta  * pts(ji,jj,jk,jp_sal) - rn_alpha * pts(ji,jj,jk,jp_tem) ) * tmask(ji,jj,jk)
225               END DO
226            END DO
227         END DO
228#else
229         DO jk = 1, jpkm1
230            prd(:,:,jk) = ( rn_beta  * pts(:,:,jk,jp_sal) - rn_alpha * pts(:,:,jk,jp_tem) ) * tmask(:,:,jk)
231         END DO
232#endif
233         !
234      END SELECT
235      !
236      IF(ln_ctl)   CALL prt_ctl( tab3d_1=prd, clinfo1=' eos  : ', ovlap=1, kdim=jpk )
237      !
238      IF( wrk_not_released(3, 1) )   CALL ctl_stop('eos_insitu: failed to release workspace array')
239      !
240
241!! * Reset control of array index permutation
242!FTRANS CLEAR
243#  include "dom_oce_ftrans.h90"
244#  include "zdfddm_ftrans.h90"
245
246   END SUBROUTINE eos_insitu
247
248
249   SUBROUTINE eos_insitu_pot( pts, prd, prhop )
250      !!----------------------------------------------------------------------
251      !!                  ***  ROUTINE eos_insitu_pot  ***
252      !!           
253      !! ** Purpose :   Compute the in situ density (ratio rho/rau0) and the
254      !!      potential volumic mass (Kg/m3) from potential temperature and
255      !!      salinity fields using an equation of state defined through the
256      !!     namelist parameter nn_eos.
257      !!
258      !! ** Method  :
259      !!      nn_eos = 0 : Jackett and McDougall (1994) equation of state.
260      !!         the in situ density is computed directly as a function of
261      !!         potential temperature relative to the surface (the opa t
262      !!         variable), salt and pressure (assuming no pressure variation
263      !!         along geopotential surfaces, i.e. the pressure p in decibars
264      !!         is approximated by the depth in meters.
265      !!              prd(t,s,p) = ( rho(t,s,p) - rau0 ) / rau0
266      !!              rhop(t,s)  = rho(t,s,0)
267      !!         with pressure                      p        decibars
268      !!              potential temperature         t        deg celsius
269      !!              salinity                      s        psu
270      !!              reference volumic mass        rau0     kg/m**3
271      !!              in situ volumic mass          rho      kg/m**3
272      !!              in situ density anomalie      prd      no units
273      !!
274      !!         Check value: rho = 1060.93298 kg/m**3 for p=10000 dbar,
275      !!          t = 40 deg celcius, s=40 psu
276      !!
277      !!      nn_eos = 1 : linear equation of state function of temperature only
278      !!              prd(t) = ( rho(t) - rau0 ) / rau0 = 0.028 - rn_alpha * t
279      !!              rhop(t,s)  = rho(t,s)
280      !!
281      !!      nn_eos = 2 : linear equation of state function of temperature and
282      !!               salinity
283      !!              prd(t,s) = ( rho(t,s) - rau0 ) / rau0
284      !!                       = rn_beta * s - rn_alpha * tn - 1.
285      !!              rhop(t,s)  = rho(t,s)
286      !!      Note that no boundary condition problem occurs in this routine
287      !!      as (tn,sn) or (ta,sa) are defined over the whole domain.
288      !!
289      !! ** Action  : - prd  , the in situ density (no units)
290      !!              - prhop, the potential volumic mass (Kg/m3)
291      !!
292      !! References :   Jackett and McDougall, J. Atmos. Ocean. Tech., 1994
293      !!                Brown and Campana, Mon. Weather Rev., 1978
294      !!----------------------------------------------------------------------
295      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
296      USE wrk_nemo, ONLY:   zws => wrk_3d_1 ! 3D workspace
297      !!
298
299!FTRANS zws :I :I :z
300!FTRANS pts :I :I :z :I
301!FTRANS prd :I :I :z
302!FTRANS prhop :I :I :z
303
304!!DCSE NEMO: This style defeats ftrans
305!     REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in   ) ::   pts    ! 1 : potential temperature  [Celcius]
306!     !                                                                ! 2 : salinity               [psu]
307!     REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(  out) ::   prd    ! in situ density            [-]
308!     REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(  out) ::   prhop  ! potential density (surface referenced)
309      REAL(wp), INTENT(in   ) ::   pts(jpi,jpj,jpkorig,jpts)   ! 1 : potential temperature  [Celcius]
310      !                                                    ! 2 : salinity               [psu]
311      REAL(wp), INTENT(  out) ::   prd(jpi,jpj,jpkorig)        ! in situ density            [-]
312      REAL(wp), INTENT(  out) ::   prhop(jpi,jpj,jpkorig)      ! potential density (surface referenced)
313      !
314      INTEGER  ::   ji, jj, jk   ! dummy loop indices
315      REAL(wp) ::   zt, zs, zh, zsr, zr1, zr2, zr3, zr4, zrhop, ze, zbw   ! local scalars
316      REAL(wp) ::   zb, zd, zc, zaw, za, zb1, za1, zkw, zk0, zrau0r       !   -      -
317      !!----------------------------------------------------------------------
318
319      IF( wrk_in_use(3, 1) ) THEN
320         CALL ctl_stop('eos_insitu_pot: requested workspace array unavailable')   ;   RETURN
321      ENDIF
322
323      SELECT CASE ( nn_eos )
324      !
325      CASE( 0 )                !==  Jackett and McDougall (1994) formulation  ==!
326         zrau0r = 1.e0 / rau0
327!CDIR NOVERRCHK
328#if defined key_z_first
329         DO jj = 1, jpj
330            DO ji = 1, jpi
331               DO jk = 1, mbkmax(ji,jj)
332                  zws(ji,jj,jk) = SQRT( ABS( pts(ji,jj,jk,jp_sal) ) )
333               END DO
334            END DO
335         END DO
336#else
337         zws(:,:,:) = SQRT( ABS( pts(:,:,:,jp_sal) ) )
338#endif
339         
340#if defined key_z_first
341         DO jj = 1, jpj
342            DO ji = 1, jpi
343               DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1
344#else
345         DO jk = 1, jpkm1
346            DO jj = 1, jpj
347               DO ji = 1, jpi
348#endif
349                  zt = pts   (ji,jj,jk,jp_tem)
350                  zs = pts   (ji,jj,jk,jp_sal)
351                  zh = fsdept(ji,jj,jk)        ! depth
352                  zsr= zws   (ji,jj,jk)        ! square root salinity
353                  !
354                  ! compute volumic mass pure water at atm pressure
355                  zr1= ( ( ( ( 6.536332e-9_wp*zt-1.120083e-6_wp )*zt+1.001685e-4_wp )*zt   &
356                     &                          -9.095290e-3_wp )*zt+6.793952e-2_wp )*zt+999.842594_wp
357                  ! seawater volumic mass atm pressure
358                  zr2= ( ( ( 5.3875e-9_wp*zt-8.2467e-7_wp ) *zt+7.6438e-5_wp ) *zt   &
359                     &                                         -4.0899e-3_wp ) *zt+0.824493_wp
360                  zr3= ( -1.6546e-6_wp*zt+1.0227e-4_wp )    *zt-5.72466e-3_wp
361                  zr4= 4.8314e-4_wp
362                  !
363                  ! potential volumic mass (reference to the surface)
364                  zrhop= ( zr4*zs + zr3*zsr + zr2 ) *zs + zr1
365                  !
366                  ! save potential volumic mass
367                  prhop(ji,jj,jk) = zrhop * tmask(ji,jj,jk)
368                  !
369                  ! add the compression terms
370                  ze = ( -3.508914e-8_wp*zt-1.248266e-8_wp ) *zt-2.595994e-6_wp
371                  zbw= (  1.296821e-6_wp*zt-5.782165e-9_wp ) *zt+1.045941e-4_wp
372                  zb = zbw + ze * zs
373                  !
374                  zd = -2.042967e-2_wp
375                  zc =   (-7.267926e-5_wp*zt+2.598241e-3_wp ) *zt+0.1571896_wp
376                  zaw= ( ( 5.939910e-6_wp*zt+2.512549e-3_wp ) *zt-0.1028859_wp ) *zt - 4.721788_wp
377                  za = ( zd*zsr + zc ) *zs + zaw
378                  !
379                  zb1=   (  -0.1909078_wp  *zt+7.390729_wp    ) *zt-55.87545_wp
380                  za1= ( (   2.326469e-3_wp*zt+1.553190_wp    ) *zt-65.00517_wp ) *zt + 1044.077_wp
381                  zkw= ( ( (-1.361629e-4_wp*zt-1.852732e-2_wp ) *zt-30.41638_wp ) *zt + 2098.925_wp ) *zt+190925.6_wp
382                  zk0= ( zb1*zsr + za1 )*zs + zkw
383                  !
384                  ! masked in situ density anomaly
385                  prd(ji,jj,jk) = (  zrhop / (  1.0_wp - zh / ( zk0 - zh * ( za - zh * zb ) )  )    &
386                     &             - rau0  ) * zrau0r * tmask(ji,jj,jk)
387               END DO
388            END DO
389         END DO
390         !
391      CASE( 1 )                !==  Linear formulation = F( temperature )  ==!
392#if defined key_z_first
393         DO jj = 1, jpj
394            DO ji = 1, jpi
395               DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1
396                  prd  (ji,jj,jk) = ( 0.0285_wp - rn_alpha * pts(ji,jj,jk,jp_tem) )        * tmask(ji,jj,jk)
397                  prhop(ji,jj,jk) = ( 1.e0_wp   +            prd(ji,jj,jk)        ) * rau0 * tmask(ji,jj,jk)
398               END DO
399            END DO
400         END DO
401#else
402         DO jk = 1, jpkm1
403            prd  (:,:,jk) = ( 0.0285_wp - rn_alpha * pts(:,:,jk,jp_tem) )        * tmask(:,:,jk)
404            prhop(:,:,jk) = ( 1.e0_wp   +            prd (:,:,jk)       ) * rau0 * tmask(:,:,jk)
405         END DO
406#endif
407         !
408      CASE( 2 )                !==  Linear formulation = F( temperature , salinity )  ==!
409#if defined key_z_first
410         DO jj = 1, jpj
411            DO ji = 1, jpi
412               DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1
413                  prd  (ji,jj,jk) = ( rn_beta  * pts(ji,jj,jk,jp_sal) - rn_alpha * pts(ji,jj,jk,jp_tem) )        * tmask(ji,jj,jk)
414                  prhop(ji,jj,jk) = ( 1.e0_wp  + prd(ji,jj,jk)                                          ) * rau0 * tmask(ji,jj,jk)
415               END DO
416            END DO
417         END DO
418#else
419         DO jk = 1, jpkm1
420            prd  (:,:,jk) = ( rn_beta  * pts(:,:,jk,jp_sal) - rn_alpha * pts(:,:,jk,jp_tem) )        * tmask(:,:,jk)
421            prhop(:,:,jk) = ( 1.e0_wp  + prd (:,:,jk) )                                       * rau0 * tmask(:,:,jk)
422         END DO
423#endif
424         !
425      END SELECT
426      !
427      IF(ln_ctl)   CALL prt_ctl( tab3d_1=prd, clinfo1=' eos-p: ', tab3d_2=prhop, clinfo2=' pot : ', ovlap=1, kdim=jpk )
428      !
429      IF( wrk_not_released(3, 1) )   CALL ctl_stop('eos_insitu_pot: failed to release workspace array')
430      !
431
432!! * Reset control of array index permutation
433!FTRANS CLEAR
434#  include "dom_oce_ftrans.h90"
435#  include "zdfddm_ftrans.h90"
436
437   END SUBROUTINE eos_insitu_pot
438
439
440   SUBROUTINE eos_insitu_2d( pts, pdep, prd )
441      !!----------------------------------------------------------------------
442      !!                  ***  ROUTINE eos_insitu_2d  ***
443      !!
444      !! ** Purpose :   Compute the in situ density (ratio rho/rau0) from
445      !!      potential temperature and salinity using an equation of state
446      !!      defined through the namelist parameter nn_eos. * 2D field case
447      !!
448      !! ** Method :
449      !!      nn_eos = 0 : Jackett and McDougall (1994) equation of state.
450      !!         the in situ density is computed directly as a function of
451      !!         potential temperature relative to the surface (the opa t
452      !!         variable), salt and pressure (assuming no pressure variation
453      !!         along geopotential surfaces, i.e. the pressure p in decibars
454      !!         is approximated by the depth in meters.
455      !!              prd(t,s,p) = ( rho(t,s,p) - rau0 ) / rau0
456      !!         with pressure                      p        decibars
457      !!              potential temperature         t        deg celsius
458      !!              salinity                      s        psu
459      !!              reference volumic mass        rau0     kg/m**3
460      !!              in situ volumic mass          rho      kg/m**3
461      !!              in situ density anomalie      prd      no units
462      !!         Check value: rho = 1060.93298 kg/m**3 for p=10000 dbar,
463      !!          t = 40 deg celcius, s=40 psu
464      !!      nn_eos = 1 : linear equation of state function of temperature only
465      !!              prd(t) = 0.0285 - rn_alpha * t
466      !!      nn_eos = 2 : linear equation of state function of temperature and
467      !!               salinity
468      !!              prd(t,s) = rn_beta * s - rn_alpha * tn - 1.
469      !!      Note that no boundary condition problem occurs in this routine
470      !!      as pts are defined over the whole domain.
471      !!
472      !! ** Action  : - prd , the in situ density (no units)
473      !!
474      !! References :   Jackett and McDougall, J. Atmos. Ocean. Tech., 1994
475      !!----------------------------------------------------------------------
476      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
477      USE wrk_nemo, ONLY:   zws => wrk_2d_5 ! 2D workspace
478      !!
479      REAL(wp), DIMENSION(jpi,jpj,jpts), INTENT(in   ) ::   pts   ! 1 : potential temperature  [Celcius]
480      !                                                           ! 2 : salinity               [psu]
481      REAL(wp), DIMENSION(jpi,jpj)     , INTENT(in   ) ::   pdep  ! depth                  [m]
482      REAL(wp), DIMENSION(jpi,jpj)     , INTENT(  out) ::   prd   ! in situ density
483      !!
484      INTEGER  ::   ji, jj                    ! dummy loop indices
485      REAL(wp) ::   zt, zs, zh, zsr, zr1, zr2, zr3, zr4, zrhop, ze, zbw   ! temporary scalars
486      REAL(wp) ::   zb, zd, zc, zaw, za, zb1, za1, zkw, zk0, zmask        !    -         -
487      !!----------------------------------------------------------------------
488
489      IF( wrk_in_use(2, 5) ) THEN
490         CALL ctl_stop('eos_insitu_2d: requested workspace array unavailable')   ;   RETURN
491      ENDIF
492
493      prd(:,:) = 0._wp
494
495      SELECT CASE( nn_eos )
496      !
497      CASE( 0 )                !==  Jackett and McDougall (1994) formulation  ==!
498      !
499!CDIR NOVERRCHK
500         DO jj = 1, jpjm1
501!CDIR NOVERRCHK
502            DO ji = 1, fs_jpim1   ! vector opt.
503               zws(ji,jj) = SQRT( ABS( pts(ji,jj,jp_sal) ) )
504            END DO
505         END DO
506         DO jj = 1, jpjm1
507            DO ji = 1, fs_jpim1   ! vector opt.
508#if defined key_z_first
509               zmask = tmask_1(ji,jj)          ! land/sea bottom mask = surf. mask
510#else
511               zmask = tmask(ji,jj,1)          ! land/sea bottom mask = surf. mask
512#endif
513               zt    = pts  (ji,jj,jp_tem)            ! interpolated T
514               zs    = pts  (ji,jj,jp_sal)            ! interpolated S
515               zsr   = zws  (ji,jj)            ! square root of interpolated S
516               zh    = pdep (ji,jj)            ! depth at the partial step level
517               !
518               ! compute volumic mass pure water at atm pressure
519               zr1 = ( ( ( ( 6.536332e-9_wp*zt-1.120083e-6_wp )*zt+1.001685e-4_wp )*zt   &
520                  &                        -9.095290e-3_wp )*zt+6.793952e-2_wp )*zt+999.842594_wp
521               ! seawater volumic mass atm pressure
522               zr2 = ( ( ( 5.3875e-9_wp*zt-8.2467e-7_wp )*zt+7.6438e-5_wp ) *zt   &
523                  &                                   -4.0899e-3_wp ) *zt+0.824493_wp
524               zr3 = ( -1.6546e-6_wp*zt+1.0227e-4_wp ) *zt-5.72466e-3_wp
525               zr4 = 4.8314e-4_wp
526               !
527               ! potential volumic mass (reference to the surface)
528               zrhop= ( zr4*zs + zr3*zsr + zr2 ) *zs + zr1
529               !
530               ! add the compression terms
531               ze = ( -3.508914e-8_wp*zt-1.248266e-8_wp ) *zt-2.595994e-6_wp
532               zbw= (  1.296821e-6_wp*zt-5.782165e-9_wp ) *zt+1.045941e-4_wp
533               zb = zbw + ze * zs
534               !
535               zd =    -2.042967e-2_wp
536               zc =   (-7.267926e-5_wp*zt+2.598241e-3_wp ) *zt+0.1571896_wp
537               zaw= ( ( 5.939910e-6_wp*zt+2.512549e-3_wp ) *zt-0.1028859_wp ) *zt -4.721788_wp
538               za = ( zd*zsr + zc ) *zs + zaw
539               !
540               zb1=     (-0.1909078_wp  *zt+7.390729_wp      ) *zt-55.87545_wp
541               za1=   ( ( 2.326469e-3_wp*zt+1.553190_wp      ) *zt-65.00517_wp ) *zt+1044.077_wp
542               zkw= ( ( (-1.361629e-4_wp*zt-1.852732e-2_wp   ) *zt-30.41638_wp ) *zt   &
543                  &                             +2098.925_wp ) *zt+190925.6_wp
544               zk0= ( zb1*zsr + za1 )*zs + zkw
545               !
546               ! masked in situ density anomaly
547               prd(ji,jj) = ( zrhop / (  1.0_wp - zh / ( zk0 - zh * ( za - zh * zb ) )  ) - rau0 ) / rau0 * zmask
548            END DO
549         END DO
550         !
551      CASE( 1 )                !==  Linear formulation = F( temperature )  ==!
552         DO jj = 1, jpjm1
553            DO ji = 1, fs_jpim1   ! vector opt.
554#if defined key_z_first
555               prd(ji,jj) = ( 0.0285_wp - rn_alpha * pts(ji,jj,jp_tem) ) * tmask_1(ji,jj)
556#else
557               prd(ji,jj) = ( 0.0285_wp - rn_alpha * pts(ji,jj,jp_tem) ) * tmask(ji,jj,1)
558#endif
559            END DO
560         END DO
561         !
562      CASE( 2 )                !==  Linear formulation = F( temperature , salinity )  ==!
563         DO jj = 1, jpjm1
564            DO ji = 1, fs_jpim1   ! vector opt.
565#if defined key_z_first
566               prd(ji,jj) = ( rn_beta * pts(ji,jj,jp_sal) - rn_alpha * pts(ji,jj,jp_tem) ) * tmask_1(ji,jj) 
567#else
568               prd(ji,jj) = ( rn_beta * pts(ji,jj,jp_sal) - rn_alpha * pts(ji,jj,jp_tem) ) * tmask(ji,jj,1) 
569#endif
570            END DO
571         END DO
572         !
573      END SELECT
574
575      IF(ln_ctl)   CALL prt_ctl( tab2d_1=prd, clinfo1=' eos2d: ' )
576      !
577      IF( wrk_not_released(2, 5) )   CALL ctl_stop('eos_insitu_2d: failed to release workspace array')
578      !
579   END SUBROUTINE eos_insitu_2d
580
581
582   SUBROUTINE eos_bn2( pts, pn2 )
583      !!----------------------------------------------------------------------
584      !!                  ***  ROUTINE eos_bn2  ***
585      !!
586      !! ** Purpose :   Compute the local Brunt-Vaisala frequency at the time-
587      !!      step of the input arguments
588      !!     
589      !! ** Method :
590      !!       * nn_eos = 0  : UNESCO sea water properties
591      !!         The brunt-vaisala frequency is computed using the polynomial
592      !!      polynomial expression of McDougall (1987):
593      !!            N^2 = grav * beta * ( alpha/beta*dk[ t ] - dk[ s ] )/e3w
594      !!      If lk_zdfddm=T, the heat/salt buoyancy flux ratio Rrau is
595      !!      computed and used in zdfddm module :
596      !!              Rrau = alpha/beta * ( dk[ t ] / dk[ s ] )
597      !!       * nn_eos = 1  : linear equation of state (temperature only)
598      !!            N^2 = grav * rn_alpha * dk[ t ]/e3w
599      !!       * nn_eos = 2  : linear equation of state (temperature & salinity)
600      !!            N^2 = grav * (rn_alpha * dk[ t ] - rn_beta * dk[ s ] ) / e3w
601      !!      The use of potential density to compute N^2 introduces e r r o r
602      !!      in the sign of N^2 at great depths. We recommand the use of
603      !!      nn_eos = 0, except for academical studies.
604      !!        Macro-tasked on horizontal slab (jk-loop)
605      !!      N.B. N^2 is set to zero at the first level (JK=1) in inidtr
606      !!      and is never used at this level.
607      !!
608      !! ** Action  : - pn2 : the brunt-vaisala frequency
609      !!
610      !! References :   McDougall, J. Phys. Oceanogr., 17, 1950-1964, 1987.
611      !!----------------------------------------------------------------------
612
613!FTRANS pts :I :I :z :I
614!FTRANS pn2 :I :I :z
615
616!!DCSE_NEMO: This style defeats ftrans
617!     REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in   ) ::   pts   ! 1 : potential temperature  [Celcius]
618!     !                                                               ! 2 : salinity               [psu]
619!     REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(  out) ::   pn2   ! Brunt-Vaisala frequency    [s-1]
620
621      REAL(wp), INTENT(in   ) ::   pts(jpi,jpj,jpkorig,jpts)   ! 1 : potential temperature  [Celcius]
622      !                                                    ! 2 : salinity               [psu]
623      REAL(wp), INTENT(  out) ::   pn2(jpi,jpj,jpkorig)        ! Brunt-Vaisala frequency    [s-1]
624      !!
625      INTEGER  ::   ji, jj, jk   ! dummy loop indices
626      REAL(wp) ::   zgde3w, zt, zs, zh, zalbet, zbeta   ! local scalars
627#if defined key_zdfddm
628      REAL(wp) ::   zds   ! local scalars
629#endif
630      !!----------------------------------------------------------------------
631
632      ! pn2 : interior points only (2=< jk =< jpkm1 )
633      ! --------------------------
634      !
635      SELECT CASE( nn_eos )
636      !
637      CASE( 0 )                !==  Jackett and McDougall (1994) formulation  ==!
638#if defined key_z_first
639         DO jj = 1, jpj
640            DO ji = 1, jpi
641               DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1
642#else
643         DO jk = 2, jpkm1
644            DO jj = 1, jpj
645               DO ji = 1, jpi
646#endif
647                  zgde3w = grav / fse3w(ji,jj,jk)
648                  zt = 0.5_wp * ( pts(ji,jj,jk,jp_tem) + pts(ji,jj,jk-1,jp_tem) )           ! potential temperature at w-pt
649                  zs = 0.5_wp * ( pts(ji,jj,jk,jp_sal) + pts(ji,jj,jk-1,jp_sal) ) - 35.0_wp ! salinity anomaly (s-35) at w-pt
650                  zh = fsdepw(ji,jj,jk)                                                ! depth in meters  at w-point
651                  !
652                  zalbet = ( ( ( - 0.255019e-07_wp * zt + 0.298357e-05_wp ) * zt   &   ! ratio alpha/beta
653                     &                                  - 0.203814e-03_wp ) * zt   &
654                     &                                  + 0.170907e-01_wp ) * zt   &
655                     &   +         0.665157e-01_wp                                 &
656                     &   +     ( - 0.678662e-05_wp * zs                            &
657                     &           - 0.846960e-04_wp * zt + 0.378110e-02_wp ) * zs   &
658                     &   +   ( ( - 0.302285e-13_wp * zh                            &
659                     &           - 0.251520e-11_wp * zs                            &
660                     &           + 0.512857e-12_wp * zt * zt              ) * zh   &
661                     &           - 0.164759e-06_wp * zs                            &
662                     &        +(   0.791325e-08_wp * zt - 0.933746e-06_wp ) * zt   &
663                     &                                  + 0.380374e-04_wp ) * zh
664                     !
665                  zbeta  = ( ( -0.415613e-09_wp * zt + 0.555579e-07_wp ) * zt      &   ! beta
666                     &                               - 0.301985e-05_wp ) * zt      &
667                     &   +       0.785567e-03_wp                                   &
668                     &   + (     0.515032e-08_wp * zs                              &
669                     &         + 0.788212e-08_wp * zt - 0.356603e-06_wp ) * zs     &
670                     &   + ( (   0.121551e-17_wp * zh                              &
671                     &         - 0.602281e-15_wp * zs                              &
672                     &         - 0.175379e-14_wp * zt + 0.176621e-12_wp ) * zh     &
673                     &                                + 0.408195e-10_wp   * zs     &
674                     &     + ( - 0.213127e-11_wp * zt + 0.192867e-09_wp ) * zt     &
675                     &                                - 0.121555e-07_wp ) * zh
676                     !
677                  pn2(ji,jj,jk) = zgde3w * zbeta * tmask(ji,jj,jk)           &   ! N^2
678                     &          * ( zalbet * ( pts(ji,jj,jk-1,jp_tem) - pts(ji,jj,jk,jp_tem) )   &
679                     &                     - ( pts(ji,jj,jk-1,jp_sal) - pts(ji,jj,jk,jp_sal) ) )
680#if defined key_zdfddm
681                  !                                                         !!bug **** caution a traiter zds=dk[S]= 0 !!!!
682                  zds = ( pts(ji,jj,jk-1,jp_sal) - pts(ji,jj,jk,jp_sal) )                    ! Rrau = (alpha / beta) (dk[t] / dk[s])
683                  IF ( ABS( zds) <= 1.e-20_wp ) zds = 1.e-20_wp
684                  rrau(ji,jj,jk) = zalbet * ( pts(ji,jj,jk-1,jp_tem) - pts(ji,jj,jk,jp_tem) ) / zds
685#endif
686               END DO
687            END DO
688         END DO
689         !
690      CASE( 1 )                !==  Linear formulation = F( temperature )  ==!
691#if defined key_z_first
692         DO jj = 1, jpj
693            DO ji = 1, jpi
694               DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1
695                  pn2(ji,jj,jk) = grav * rn_alpha * ( pts(ji,jj,jk-1,jp_tem) - pts(ji,jj,jk,jp_tem) ) &
696                     &                 / fse3w(ji,jj,jk) * tmask(ji,jj,jk)
697               END DO
698            END DO
699         END DO
700#else
701         DO jk = 2, jpkm1
702            pn2(:,:,jk) = grav * rn_alpha * ( pts(:,:,jk-1,jp_tem) - pts(:,:,jk,jp_tem) ) / fse3w(:,:,jk) * tmask(:,:,jk)
703         END DO
704#endif
705         !
706      CASE( 2 )                !==  Linear formulation = F( temperature , salinity )  ==!
707#if defined key_z_first
708         DO jj = 1, jpj
709            DO ji = 1, jpi
710               DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1
711                  pn2(ji,jj,jk) = grav * (  rn_alpha * ( pts(ji,jj,jk-1,jp_tem) - pts(ji,jj,jk,jp_tem) )      &
712                     &                    - rn_beta  * ( pts(ji,jj,jk-1,jp_sal) - pts(ji,jj,jk,jp_sal) )  )   &
713                     &                 / fse3w(ji,jj,jk) * tmask(ji,jj,jk)
714               END DO
715            END DO
716         END DO 
717#else
718         DO jk = 2, jpkm1
719            pn2(:,:,jk) = grav * (  rn_alpha * ( pts(:,:,jk-1,jp_tem) - pts(:,:,jk,jp_tem) )      &
720               &                  - rn_beta  * ( pts(:,:,jk-1,jp_sal) - pts(:,:,jk,jp_sal) )  )   &
721               &               / fse3w(:,:,jk) * tmask(:,:,jk)
722         END DO 
723#endif
724#if defined key_zdfddm
725#if defined key_z_first
726         DO jj = 1, jpj                                   ! Rrau = (alpha / beta) (dk[t] / dk[s])
727            DO ji = 1, jpi
728               DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1
729#else
730         DO jk = 2, jpkm1                                 ! Rrau = (alpha / beta) (dk[t] / dk[s])
731            DO jj = 1, jpj
732               DO ji = 1, jpi
733#endif
734                  zds = ( pts(ji,jj,jk-1,jp_sal) - pts(ji,jj,jk,jp_sal) ) 
735                  IF ( ABS( zds ) <= 1.e-20_wp ) zds = 1.e-20_wp
736                  rrau(ji,jj,jk) = ralpbet * ( pts(ji,jj,jk-1,jp_tem) - pts(ji,jj,jk,jp_tem) ) / zds
737               END DO
738            END DO
739         END DO
740#endif
741      END SELECT
742
743      IF(ln_ctl)   CALL prt_ctl( tab3d_1=pn2, clinfo1=' bn2  : ', ovlap=1, kdim=jpk )
744#if defined key_zdfddm
745      IF(ln_ctl)   CALL prt_ctl( tab3d_1=rrau, clinfo1=' rrau : ', ovlap=1, kdim=jpk )
746#endif
747      !
748
749!! * Reset control of array index permutation
750!FTRANS CLEAR
751#  include "dom_oce_ftrans.h90"
752#  include "zdfddm_ftrans.h90"
753
754   END SUBROUTINE eos_bn2
755
756
757   SUBROUTINE eos_alpbet( pts, palph, pbeta )
758      !!----------------------------------------------------------------------
759      !!                 ***  ROUTINE ldf_slp_grif  ***
760      !!
761      !! ** Purpose :   Calculates the thermal and haline expansion coefficients at T-points.
762      !!
763      !! ** Method  :   calculates alpha and beta at T-points
764      !!       * nn_eos = 0  : UNESCO sea water properties
765      !!         The brunt-vaisala frequency is computed using the polynomial
766      !!      polynomial expression of McDougall (1987):
767      !!            N^2 = grav * beta * ( alpha/beta*dk[ t ] - dk[ s ] )/e3w
768      !!      If lk_zdfddm=T, the heat/salt buoyancy flux ratio Rrau is
769      !!      computed and used in zdfddm module :
770      !!              Rrau = alpha/beta * ( dk[ t ] / dk[ s ] )
771      !!       * nn_eos = 1  : linear equation of state (temperature only)
772      !!            N^2 = grav * rn_alpha * dk[ t ]/e3w
773      !!       * nn_eos = 2  : linear equation of state (temperature & salinity)
774      !!            N^2 = grav * (rn_alpha * dk[ t ] - rn_beta * dk[ s ] ) / e3w
775      !!       * nn_eos = 3  : Jackett JAOT 2003 ???
776      !!
777      !! ** Action  : - palph, pbeta : thermal and haline expansion coeff. at T-point
778      !!----------------------------------------------------------------------
779
780!FTRANS pts :I :I :z :I
781!FTRANS palph :I :I :z
782!FTRANS pbeta :I :I :z
783!!DCSE_NEMO: This style defeats ftrans
784!     REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in   ) ::   pts            ! pot. temperature & salinity
785!     REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(  out) ::   palph, pbeta   ! thermal & haline expansion coeff.
786      REAL(wp), INTENT(in   ) ::   pts(jpi,jpj,jpkorig,jpts)            ! pot. temperature & salinity
787      REAL(wp), INTENT(  out) ::   palph(jpi,jpj,jpkorig)               ! thermal expansion coeff.
788      REAL(wp), INTENT(  out) ::   pbeta(jpi,jpj,jpkorig)               ! haline  expansion coeff.
789      !
790      INTEGER  ::   ji, jj, jk   ! dummy loop indices
791      REAL(wp) ::   zt, zs, zh   ! local scalars
792      !!----------------------------------------------------------------------
793      !
794      SELECT CASE ( nn_eos )
795      !
796      CASE ( 0 )               ! Jackett and McDougall (1994) formulation
797#if defined key_z_first
798         DO jj = 1, jpj
799            DO ji = 1, jpi
800               DO jk = 1, mbkmax(ji,jj) ! jpk
801#else
802         DO jk = 1, jpk
803            DO jj = 1, jpj
804               DO ji = 1, jpi
805#endif
806                  zt = pts(ji,jj,jk,jp_tem)           ! potential temperature
807                  zs = pts(ji,jj,jk,jp_sal) - 35._wp  ! salinity anomaly (s-35)
808                  zh = fsdept(ji,jj,jk)              ! depth in meters
809                  !
810                  pbeta(ji,jj,jk) = ( ( -0.415613e-09_wp * zt + 0.555579e-07_wp ) * zt   &
811                     &                                        - 0.301985e-05_wp ) * zt   &
812                     &           + 0.785567e-03_wp                                       &
813                     &           + (     0.515032e-08_wp * zs                            &
814                     &                 + 0.788212e-08_wp * zt - 0.356603e-06_wp ) * zs   &
815                     &           + ( (   0.121551e-17_wp * zh                            &
816                     &                 - 0.602281e-15_wp * zs                            &
817                     &                 - 0.175379e-14_wp * zt + 0.176621e-12_wp ) * zh   &
818                     &                                        + 0.408195e-10_wp   * zs   &
819                     &             + ( - 0.213127e-11_wp * zt + 0.192867e-09_wp ) * zt   &
820                     &                                        - 0.121555e-07_wp ) * zh
821                     !
822                  palph(ji,jj,jk) = - pbeta(ji,jj,jk) *                             &
823                      &     ((( ( - 0.255019e-07_wp * zt + 0.298357e-05_wp ) * zt   &
824                      &                                  - 0.203814e-03_wp ) * zt   &
825                      &                                  + 0.170907e-01_wp ) * zt   &
826                      &   + 0.665157e-01_wp                                         &
827                      &   +     ( - 0.678662e-05_wp * zs                            &
828                      &           - 0.846960e-04_wp * zt + 0.378110e-02_wp ) * zs   &
829                      &   +   ( ( - 0.302285e-13_wp * zh                            &
830                      &           - 0.251520e-11_wp * zs                            &
831                      &           + 0.512857e-12_wp * zt * zt              ) * zh   &
832                      &           - 0.164759e-06_wp * zs                            &
833                      &        +(   0.791325e-08_wp * zt - 0.933746e-06_wp ) * zt   &
834                      &                                  + 0.380374e-04_wp ) * zh)
835               END DO
836            END DO
837         END DO
838         !
839      CASE ( 1 )
840         palph(:,:,:) = - rn_alpha
841         pbeta(:,:,:) =   0._wp
842         !
843      CASE ( 2 )
844         palph(:,:,:) = - rn_alpha
845         pbeta(:,:,:) =   rn_beta
846         !
847      CASE DEFAULT
848         IF(lwp) WRITE(numout,cform_err)
849         IF(lwp) WRITE(numout,*) '          bad flag value for nn_eos = ', nn_eos
850         nstop = nstop + 1
851         !
852      END SELECT
853      !
854
855!! * Reset control of array index permutation
856!FTRANS CLEAR
857#  include "dom_oce_ftrans.h90"
858#  include "zdfddm_ftrans.h90"
859
860   END SUBROUTINE eos_alpbet
861
862
863   FUNCTION tfreez( psal ) RESULT( ptf )
864      !!----------------------------------------------------------------------
865      !!                 ***  ROUTINE eos_init  ***
866      !!
867      !! ** Purpose :   Compute the sea surface freezing temperature [Celcius]
868      !!
869      !! ** Method  :   UNESCO freezing point at the surface (pressure = 0???)
870      !!       freezing point [Celcius]=(-.0575+1.710523e-3*sqrt(abs(s))-2.154996e-4*s)*s-7.53e-4*p
871      !!       checkvalue: tf= -2.588567 Celsius for s=40.0psu, p=500. decibars
872      !!
873      !! Reference  :   UNESCO tech. papers in the marine science no. 28. 1978
874      !!----------------------------------------------------------------------
875      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ) ::   psal   ! salinity             [psu]
876      ! Leave result array automatic rather than making explicitly allocated
877      REAL(wp), DIMENSION(jpi,jpj)                ::   ptf    ! freezing temperature [Celcius]
878      !!----------------------------------------------------------------------
879      !
880      ptf(:,:) = ( - 0.0575_wp + 1.710523e-3_wp * SQRT( psal(:,:) )   &
881         &                     - 2.154996e-4_wp *       psal(:,:)   ) * psal(:,:)
882      !
883   END FUNCTION tfreez
884
885
886   SUBROUTINE eos_init
887      !!----------------------------------------------------------------------
888      !!                 ***  ROUTINE eos_init  ***
889      !!
890      !! ** Purpose :   initializations for the equation of state
891      !!
892      !! ** Method  :   Read the namelist nameos and control the parameters
893      !!----------------------------------------------------------------------
894      NAMELIST/nameos/ nn_eos, rn_alpha, rn_beta
895      !!----------------------------------------------------------------------
896      !
897      REWIND( numnam )            ! Read Namelist nameos : equation of state
898      READ  ( numnam, nameos )
899      !
900      IF(lwp) THEN                ! Control print
901         WRITE(numout,*)
902         WRITE(numout,*) 'eos_init : equation of state'
903         WRITE(numout,*) '~~~~~~~~'
904         WRITE(numout,*) '          Namelist nameos : set eos parameters'
905         WRITE(numout,*) '             flag for eq. of state and N^2  nn_eos   = ', nn_eos
906         WRITE(numout,*) '             thermal exp. coef. (linear)    rn_alpha = ', rn_alpha
907         WRITE(numout,*) '             saline  exp. coef. (linear)    rn_beta  = ', rn_beta
908      ENDIF
909      !
910      SELECT CASE( nn_eos )         ! check option
911      !
912      CASE( 0 )                        !==  Jackett and McDougall (1994) formulation  ==!
913         IF(lwp) WRITE(numout,*)
914         IF(lwp) WRITE(numout,*) '          use of Jackett & McDougall (1994) equation of state and'
915         IF(lwp) WRITE(numout,*) '                 McDougall (1987) Brunt-Vaisala frequency'
916         !
917      CASE( 1 )                        !==  Linear formulation = F( temperature )  ==!
918         IF(lwp) WRITE(numout,*)
919         IF(lwp) WRITE(numout,*) '          use of linear eos rho(T) = rau0 * ( 1.0285 - rn_alpha * T )'
920         IF( lk_zdfddm ) CALL ctl_stop( '          double diffusive mixing parameterization requires',   &
921              &                         ' that T and S are used as state variables' )
922         !
923      CASE( 2 )                        !==  Linear formulation = F( temperature , salinity )  ==!
924         ralpbet = rn_alpha / rn_beta
925         IF(lwp) WRITE(numout,*)
926         IF(lwp) WRITE(numout,*) '          use of linear eos rho(T,S) = rau0 * ( rn_beta * S - rn_alpha * T )'
927         !
928      CASE DEFAULT                     !==  ERROR in nn_eos  ==!
929         WRITE(ctmp1,*) '          bad flag value for nn_eos = ', nn_eos
930         CALL ctl_stop( ctmp1 )
931         !
932      END SELECT
933      !
934   END SUBROUTINE eos_init
935
936   !!======================================================================
937END MODULE eosbn2 
Note: See TracBrowser for help on using the repository browser.