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.
cool_skin.F90 in NEMO/branches/UKMO/dev_r9885_GO6_mixing/src/OCE/DIU – NEMO

source: NEMO/branches/UKMO/dev_r9885_GO6_mixing/src/OCE/DIU/cool_skin.F90 @ 9895

Last change on this file since 9895 was 9895, checked in by davestorkey, 6 years ago

UKMO/dev_r9885_GO6_mixing branch: clear SVN keywords.

File size: 6.5 KB
Line 
1MODULE cool_skin
2   !!======================================================================
3   !!                    ***  MODULE  cool_skin  ***
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
41   PUBLIC diurnal_sst_coolskin_step, diurnal_sst_coolskin_init
42
43      !! * Substitutions
44#  include "vectopt_loop_substitute.h90"
45   !!----------------------------------------------------------------------
46   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
47   !! $Id$
48   !! Software governed by the CeCILL licence     (./LICENSE)
49   !!----------------------------------------------------------------------   
50   CONTAINS
51   
52   SUBROUTINE diurnal_sst_coolskin_init
53      !!----------------------------------------------------------------------
54      !! *** ROUTINE diurnal_sst_coolskin_init ***
55      !!
56      !! ** Purpose :   initialise the cool skin model
57      !!
58      !! ** Method :
59      !!
60      !! ** Reference :
61      !!
62      !!----------------------------------------------------------------------
63      ALLOCATE( x_csdsst(jpi,jpj), x_csthick(jpi,jpj) )
64      x_csdsst = 0.
65      x_csthick = 0.
66      !
67   END SUBROUTINE diurnal_sst_coolskin_init
68
69
70   SUBROUTINE diurnal_sst_coolskin_step(psqflux, pstauflux, psrho, rdt)
71      !!----------------------------------------------------------------------
72      !! *** ROUTINE diurnal_sst_takaya_step ***
73      !!
74      !! ** Purpose :   Time-step the Artale cool skin model
75      !!
76      !! ** Method :
77      !!
78      !! ** Reference :
79      !!----------------------------------------------------------------------
80      ! Dummy variables
81      REAL(wp), INTENT(IN), DIMENSION(jpi,jpj) :: psqflux     ! Heat (non-solar)(Watts)
82      REAL(wp), INTENT(IN), DIMENSION(jpi,jpj) :: pstauflux   ! Wind stress (kg/ m s^2)
83      REAL(wp), INTENT(IN), DIMENSION(jpi,jpj) :: psrho       ! Water density (kg/m^3)
84      REAL(wp), INTENT(IN) :: rdt                             ! Time-step
85     
86      ! Local variables
87      REAL(wp), DIMENSION(jpi,jpj) :: z_fv                    ! Friction velocity     
88      REAL(wp), DIMENSION(jpi,jpj) :: z_gamma                 ! Dimensionless function of wind speed
89      REAL(wp), DIMENSION(jpi,jpj) :: z_lamda                 ! Sauders (dimensionless) proportionality constant
90      REAL(wp), DIMENSION(jpi,jpj) :: z_wspd                  ! Wind speed (m/s)
91      REAL(wp) :: z_ztx                                       ! Temporary u wind stress
92      REAL(wp) :: z_zty                                       ! Temporary v wind stress
93      REAL(wp) :: z_zmod                                      ! Temporary total wind stress
94     
95      INTEGER :: ji,jj
96      !!----------------------------------------------------------------------
97      !
98      IF( .NOT. ln_blk )   CALL ctl_stop("cool_skin.f90: diurnal flux processing only implemented for bulk forcing")
99      !
100      DO jj = 1,jpj
101         DO ji = 1,jpi
102            !
103            ! Calcualte wind speed from wind stress and friction velocity
104            IF( tmask(ji,jj,1) == 1. .AND. pstauflux(ji,jj) /= 0 .AND. psrho(ji,jj) /=0 ) THEN
105               z_fv(ji,jj) = SQRT( pstauflux(ji,jj) / psrho(ji,jj) )
106               z_wspd(ji,jj) = SQRT( pstauflux(ji,jj) / ( pp_cda * pp_rhoa ) )
107            ELSE
108               z_fv(ji,jj) = 0.
109               z_wspd(ji,jj) = 0.     
110            ENDIF
111            !
112            ! Calculate gamma function which is dependent upon wind speed
113            IF( tmask(ji,jj,1) == 1. ) THEN
114               IF( ( z_wspd(ji,jj) <= 7.5 ) ) z_gamma(ji,jj) = ( 0.2 * z_wspd(ji,jj) ) + 0.5
115               IF( ( z_wspd(ji,jj) > 7.5 ) .AND. ( z_wspd(ji,jj) < 10. ) ) z_gamma(ji,jj) = ( 1.6 * z_wspd(ji,jj) ) - 10.
116               IF( ( z_wspd(ji,jj) >= 10. ) ) z_gamma(ji,jj) = 6.
117            ENDIF
118            !
119            ! Calculate lamda function
120            IF( tmask(ji,jj,1) == 1. .AND. z_fv(ji,jj) /= 0 ) THEN
121               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 )
122            ELSE
123               z_lamda(ji,jj) = 0.
124            ENDIF
125            !
126            ! Calculate the cool skin thickness - only when heat flux is out of the ocean
127            IF( tmask(ji,jj,1) == 1. .AND. z_fv(ji,jj) /= 0 .AND. psqflux(ji,jj) < 0 ) THEN
128               x_csthick(ji,jj) = ( z_lamda(ji,jj) * pp_v ) / z_fv(ji,jj)
129            ELSE
130               x_csthick(ji,jj) = 0.
131            ENDIF
132            !
133            ! Calculate the cool skin correction - only when the heat flux is out of the ocean
134            IF( tmask(ji,jj,1) == 1. .AND. x_csthick(ji,jj) /= 0. .AND. psqflux(ji,jj) < 0. ) THEN
135               x_csdsst(ji,jj) = ( psqflux(ji,jj) * x_csthick(ji,jj) ) / pp_k
136             ELSE
137               x_csdsst(ji,jj) = 0.
138            ENDIF
139            !
140         END DO
141      END DO
142      !
143   END SUBROUTINE diurnal_sst_coolskin_step
144
145   !!======================================================================
146END MODULE cool_skin
Note: See TracBrowser for help on using the repository browser.