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 branches/UKMO/AMM15_v3_6_STABLE_package/NEMOGCM/NEMO/OPA_SRC/DIU – NEMO

source: branches/UKMO/AMM15_v3_6_STABLE_package/NEMOGCM/NEMO/OPA_SRC/DIU/cool_skin.F90 @ 8561

Last change on this file since 8561 was 8561, checked in by jgraham, 7 years ago

Updates for operational diagnostics:
25h mean diagnostics - bottom temperature (and insitu temp)
Operational foam diagnostics - diaopfoam and DIU routines added.

File size: 5.9 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_step  : timestep the cool skin corrections
13   !!----------------------------------------------------------------------
14   USE par_kind
15   USE phycst
16   USE dom_oce
17   USE in_out_manager
18   USE sbc_oce
19   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
20   
21   IMPLICIT NONE
22   
23   !Namelist parameters
24
25   !Parameters
26   REAL(wp), PRIVATE, PARAMETER :: pp_k = 0.596_wp          !Thermal conductivity of seawater
27   REAL(wp), PRIVATE, PARAMETER :: pp_v = 1.05e-6_wp        !Kinematic viscosity of seawater
28   REAL(wp), PRIVATE, PARAMETER :: pp_C = 86400             !seconds [see Tu and Tsuang (2005)]
29   REAL(wp), PRIVATE, PARAMETER :: pp_cw = 3993._wp         !specific heat capacity of seawater
30   REAL(wp), PRIVATE, PARAMETER :: pp_h = 10._wp            !reference depth [using 10m from Artale et al. (2002)]
31   REAL(wp), PRIVATE, PARAMETER :: pp_rhoa = 1.20421_wp     !density of air (at 20C)
32   REAL(wp), PRIVATE, PARAMETER :: pp_cda = 1.45e-3_wp      !assumed air-sea drag coefficient for calculating wind speed
33   
34   !Key variables
35   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: x_csdsst     !Cool skin delta SST
36   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: x_csthick  !Cool skin thickness
37
38   PRIVATE
39   PUBLIC diurnal_sst_coolskin_step, diurnal_sst_coolskin_init
40
41      !! * Substitutions
42#  include "domzgr_substitute.h90"
43#  include "vectopt_loop_substitute.h90"
44   
45   CONTAINS
46   
47   SUBROUTINE diurnal_sst_coolskin_init
48      !!----------------------------------------------------------------------
49      !! *** ROUTINE diurnal_sst_coolskin_init ***
50      !!
51      !! ** Purpose :   initalise the coolskin model
52      !!
53      !! ** Method :
54      !!
55      !! ** Reference :
56      !!
57      !!----------------------------------------------------------------------
58     
59      IMPLICIT NONE
60     
61      ALLOCATE( x_csdsst(jpi,jpj), x_csthick(jpi,jpj) )
62      x_csdsst = 0.
63      x_csthick = 0.
64     
65   END SUBROUTINE diurnal_sst_coolskin_init
66 
67   SUBROUTINE diurnal_sst_coolskin_step(psqflux, pstauflux, psrho, rdt)
68      !!----------------------------------------------------------------------
69      !! *** ROUTINE diurnal_sst_takaya_step ***
70      !!
71      !! ** Purpose :   Timestep the Artale cool skin model
72      !!
73      !! ** Method :
74      !!
75      !! ** Reference :
76      !!----------------------------------------------------------------------
77     
78      IMPLICIT NONE
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                             !Timestep
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 windspeed
89      REAL(wp), DIMENSION(jpi,jpj) :: z_lamda                 !Sauders (dimensionless) proportionality constant
90      REAL(wp), DIMENSION(jpi,jpj) :: z_wspd                  !Windspeed (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      DO jj = 1,jpj
98         DO ji = 1,jpi
99           
100            !Calcualte wind speed from wind stress and friction velocity
101            IF( tmask(ji,jj,1) == 1. .AND. pstauflux(ji,jj) /= 0 .AND. psrho(ji,jj) /=0 ) THEN
102               z_fv(ji,jj) = SQRT( pstauflux(ji,jj) / psrho(ji,jj) )
103               z_wspd(ji,jj) = SQRT( pstauflux(ji,jj) / ( pp_cda * pp_rhoa ) )
104            ELSE
105               z_fv(ji,jj) = 0.
106               z_wspd(ji,jj) = 0.     
107            ENDIF
108
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
118            !Calulate lamda function
119            IF( tmask(ji,jj,1) == 1. .AND. z_fv(ji,jj) /= 0 ) THEN
120               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 )
121            ELSE
122               z_lamda(ji,jj) = 0.
123            ENDIF
124
125
126
127            !Calculate the cool skin thickness - only when heat flux is out of the ocean
128            IF( tmask(ji,jj,1) == 1. .AND. z_fv(ji,jj) /= 0 .AND. psqflux(ji,jj) < 0 ) THEN
129                x_csthick(ji,jj) = ( z_lamda(ji,jj) * pp_v ) / z_fv(ji,jj)
130            ELSE
131                x_csthick(ji,jj) = 0.
132            ENDIF
133
134
135
136            !Calculate the cool skin correction - only when the heat flux is out of the ocean
137            IF( tmask(ji,jj,1) == 1. .AND. x_csthick(ji,jj) /= 0. .AND. psqflux(ji,jj) < 0. ) THEN
138               x_csdsst(ji,jj) = ( psqflux(ji,jj) * x_csthick(ji,jj) ) / pp_k
139             ELSE
140               x_csdsst(ji,jj) = 0.
141            ENDIF
142
143         ENDDO
144      ENDDO
145
146   END SUBROUTINE diurnal_sst_coolskin_step
147
148
149END MODULE cool_skin
Note: See TracBrowser for help on using the repository browser.