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.
diu_coolskin.F90 in NEMO/trunk/src/OCE/DIU – NEMO

source: NEMO/trunk/src/OCE/DIU/diu_coolskin.F90 @ 14072

Last change on this file since 14072 was 14072, checked in by laurent, 3 years ago

Merging branch "2020/dev_r13648_ASINTER-04_laurent_bulk_ice", ticket #2369

  • Property svn:keywords set to Id
File size: 6.4 KB
Line 
1MODULE diu_coolskin
2   !!======================================================================
3   !!                    ***  MODULE  diu_coolskin  ***
4   !!     Cool skin thickness and delta T correction using Artele et al. (2002)
5   !!     [see also Tu and Tsuang (2005)]
6   !!
7   !!=====================================================================
8   !! History :        !  2012-01  (P. Sykes)  Original code
9   !!----------------------------------------------------------------------
10
11   !!----------------------------------------------------------------------
12   !!   diurnal_sst_coolskin_init : initialisation of the cool skin
13   !!   diurnal_sst_coolskin_step : time-stepping  of the cool skin corrections
14   !!----------------------------------------------------------------------
15   USE par_kind
16   USE phycst
17   USE dom_oce
18   USE in_out_manager
19   USE sbc_oce
20   USE lib_mpp
21   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
22
23   IMPLICIT NONE
24   PRIVATE
25
26   ! Namelist parameters
27
28   ! Parameters
29   REAL(wp), PRIVATE, PARAMETER :: pp_k = 0.596_wp          ! Thermal conductivity of seawater
30   REAL(wp), PRIVATE, PARAMETER :: pp_v = 1.05e-6_wp        ! Kinematic viscosity of seawater
31   REAL(wp), PRIVATE, PARAMETER :: pp_C = 86400             ! seconds [see Tu and Tsuang (2005)]
32   REAL(wp), PRIVATE, PARAMETER :: pp_cw = 3993._wp         ! specific heat capacity of seawater
33   REAL(wp), PRIVATE, PARAMETER :: pp_h = 10._wp            ! reference depth [using 10m from Artale et al. (2002)]
34   REAL(wp), PRIVATE, PARAMETER :: pp_rhoa = 1.20421_wp     ! density of air (at 20C)
35   REAL(wp), PRIVATE, PARAMETER :: pp_cda = 1.45e-3_wp      ! assumed air-sea drag coefficient for calculating wind speed
36
37   ! Key variables
38   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: x_csdsst    ! Cool skin delta SST
39   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: x_csthick   ! Cool skin thickness
40   PUBLIC diurnal_sst_coolskin_step, diurnal_sst_coolskin_init
41
42      !! * Substitutions
43#  include "do_loop_substitute.h90"
44   !!----------------------------------------------------------------------
45   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
46   !! $Id$
47   !! Software governed by the CeCILL license (see ./LICENSE)
48   !!----------------------------------------------------------------------
49   CONTAINS
50
51   SUBROUTINE diurnal_sst_coolskin_init
52      !!----------------------------------------------------------------------
53      !! *** ROUTINE diurnal_sst_coolskin_init ***
54      !!
55      !! ** Purpose :   initialise the cool skin model
56      !!
57      !! ** Method :
58      !!
59      !! ** Reference :
60      !!
61      !!----------------------------------------------------------------------
62      ALLOCATE( x_csdsst(jpi,jpj), x_csthick(jpi,jpj) )
63      x_csdsst = 0.
64      x_csthick = 0.
65      !
66   END SUBROUTINE diurnal_sst_coolskin_init
67
68
69   SUBROUTINE diurnal_sst_coolskin_step(psqflux, pstauflux, psrho, pDt)
70      !!----------------------------------------------------------------------
71      !! *** ROUTINE diurnal_sst_takaya_step ***
72      !!
73      !! ** Purpose :   Time-step the Artale cool skin model
74      !!
75      !! ** Method :
76      !!
77      !! ** Reference :
78      !!----------------------------------------------------------------------
79      ! Dummy variables
80      REAL(wp), INTENT(IN), DIMENSION(jpi,jpj) :: psqflux     ! Heat (non-solar)(Watts)
81      REAL(wp), INTENT(IN), DIMENSION(jpi,jpj) :: pstauflux   ! Wind stress (kg/ m s^2)
82      REAL(wp), INTENT(IN), DIMENSION(jpi,jpj) :: psrho       ! Water density (kg/m^3)
83      REAL(wp), INTENT(IN) :: pDt                             ! Time-step
84
85      ! Local variables
86      REAL(wp), DIMENSION(jpi,jpj) :: z_fv                    ! Friction velocity
87      REAL(wp), DIMENSION(jpi,jpj) :: z_gamma                 ! Dimensionless function of wind speed
88      REAL(wp), DIMENSION(jpi,jpj) :: z_lamda                 ! Sauders (dimensionless) proportionality constant
89      REAL(wp), DIMENSION(jpi,jpj) :: z_wspd                  ! Wind speed (m/s)
90      REAL(wp) :: z_ztx                                       ! Temporary u wind stress
91      REAL(wp) :: z_zty                                       ! Temporary v wind stress
92      REAL(wp) :: z_zmod                                      ! Temporary total wind stress
93
94      INTEGER :: ji,jj
95      !!----------------------------------------------------------------------
96      !
97      IF( .NOT. (ln_blk .OR. ln_abl) )   CALL ctl_stop("diu_coolskin.f90: diurnal flux processing only implemented for bulk forcing")
98      !
99      DO_2D( 1, 1, 1, 1 )
100         !
101         ! Calcualte wind speed from wind stress and friction velocity
102         IF( tmask(ji,jj,1) == 1. .AND. pstauflux(ji,jj) /= 0 .AND. psrho(ji,jj) /=0 ) THEN
103            z_fv(ji,jj) = SQRT( pstauflux(ji,jj) / psrho(ji,jj) )
104            z_wspd(ji,jj) = SQRT( pstauflux(ji,jj) / ( pp_cda * pp_rhoa ) )
105         ELSE
106            z_fv(ji,jj) = 0.
107            z_wspd(ji,jj) = 0.
108         ENDIF
109         !
110         ! Calculate gamma function which is dependent upon wind speed
111         IF( tmask(ji,jj,1) == 1. ) THEN
112            IF( ( z_wspd(ji,jj) <= 7.5 ) ) z_gamma(ji,jj) = ( 0.2 * z_wspd(ji,jj) ) + 0.5
113            IF( ( z_wspd(ji,jj) > 7.5 ) .AND. ( z_wspd(ji,jj) < 10. ) ) z_gamma(ji,jj) = ( 1.6 * z_wspd(ji,jj) ) - 10.
114            IF( ( z_wspd(ji,jj) >= 10. ) ) z_gamma(ji,jj) = 6.
115         ENDIF
116         !
117         ! Calculate lamda function
118         IF( tmask(ji,jj,1) == 1. .AND. z_fv(ji,jj) /= 0 ) THEN
119            z_lamda(ji,jj) = ( z_fv(ji,jj) * pp_k * pp_C ) / ( z_gamma(ji,jj) * psrho(ji,jj) * pp_cw * pp_h * pp_v )
120         ELSE
121            z_lamda(ji,jj) = 0.
122         ENDIF
123         !
124         ! Calculate the cool skin thickness - only when heat flux is out of the ocean
125         IF( tmask(ji,jj,1) == 1. .AND. z_fv(ji,jj) /= 0 .AND. psqflux(ji,jj) < 0 ) THEN
126            x_csthick(ji,jj) = ( z_lamda(ji,jj) * pp_v ) / z_fv(ji,jj)
127         ELSE
128            x_csthick(ji,jj) = 0.
129         ENDIF
130         !
131         ! Calculate the cool skin correction - only when the heat flux is out of the ocean
132         IF( tmask(ji,jj,1) == 1. .AND. x_csthick(ji,jj) /= 0. .AND. psqflux(ji,jj) < 0. ) THEN
133            x_csdsst(ji,jj) = ( psqflux(ji,jj) * x_csthick(ji,jj) ) / pp_k
134          ELSE
135            x_csdsst(ji,jj) = 0.
136         ENDIF
137         !
138      END_2D
139      !
140   END SUBROUTINE diurnal_sst_coolskin_step
141
142   !!======================================================================
143END MODULE diu_coolskin
Note: See TracBrowser for help on using the repository browser.