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/2015/dev_merge_2015/NEMOGCM/NEMO/OPA_SRC/DIU – NEMO

source: branches/2015/dev_merge_2015/NEMOGCM/NEMO/OPA_SRC/DIU/cool_skin.F90 @ 6075

Last change on this file since 6075 was 6075, checked in by timgraham, 8 years ago

Lots of bug fixes to get sette to compile

File size: 6.1 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  : time-step 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 "vectopt_loop_substitute.h90"
43   
44   CONTAINS
45   
46   SUBROUTINE diurnal_sst_coolskin_init
47      !!----------------------------------------------------------------------
48      !! *** ROUTINE diurnal_sst_coolskin_init ***
49      !!
50      !! ** Purpose :   initialise the cool skin model
51      !!
52      !! ** Method :
53      !!
54      !! ** Reference :
55      !!
56      !!----------------------------------------------------------------------
57     
58      IMPLICIT NONE
59     
60      ALLOCATE( x_csdsst(jpi,jpj), x_csthick(jpi,jpj) )
61      x_csdsst = 0.
62      x_csthick = 0.
63     
64   END SUBROUTINE diurnal_sst_coolskin_init
65 
66   SUBROUTINE diurnal_sst_coolskin_step(psqflux, pstauflux, psrho, rdt)
67      !!----------------------------------------------------------------------
68      !! *** ROUTINE diurnal_sst_takaya_step ***
69      !!
70      !! ** Purpose :   Time-step the Artale cool skin model
71      !!
72      !! ** Method :
73      !!
74      !! ** Reference :
75      !!----------------------------------------------------------------------
76     
77      IMPLICIT NONE
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) :: rdt                             ! 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      IF ( .NOT. ln_blk_core ) THEN
97         CALL ctl_stop("cool_skin.f90: diurnal flux processing only implemented"//&
98         &             " for core bulk forcing")
99      ENDIF
100 
101      DO jj = 1,jpj
102         DO ji = 1,jpi
103           
104            ! Calcualte wind speed from wind stress and friction velocity
105            IF( tmask(ji,jj,1) == 1. .AND. pstauflux(ji,jj) /= 0 .AND. psrho(ji,jj) /=0 ) THEN
106               z_fv(ji,jj) = SQRT( pstauflux(ji,jj) / psrho(ji,jj) )
107               z_wspd(ji,jj) = SQRT( pstauflux(ji,jj) / ( pp_cda * pp_rhoa ) )
108            ELSE
109               z_fv(ji,jj) = 0.
110               z_wspd(ji,jj) = 0.     
111            ENDIF
112
113 
114            ! Calculate gamma function which is dependent upon wind speed
115            IF( tmask(ji,jj,1) == 1. ) THEN
116               IF( ( z_wspd(ji,jj) <= 7.5 ) ) z_gamma(ji,jj) = ( 0.2 * z_wspd(ji,jj) ) + 0.5
117               IF( ( z_wspd(ji,jj) > 7.5 ) .AND. ( z_wspd(ji,jj) < 10. ) ) z_gamma(ji,jj) = ( 1.6 * z_wspd(ji,jj) ) - 10.
118               IF( ( z_wspd(ji,jj) >= 10. ) ) z_gamma(ji,jj) = 6.
119            ENDIF
120
121
122            ! Calculate lamda function
123            IF( tmask(ji,jj,1) == 1. .AND. z_fv(ji,jj) /= 0 ) THEN
124               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 )
125            ELSE
126               z_lamda(ji,jj) = 0.
127            ENDIF
128
129
130
131            ! Calculate the cool skin thickness - only when heat flux is out of the ocean
132            IF( tmask(ji,jj,1) == 1. .AND. z_fv(ji,jj) /= 0 .AND. psqflux(ji,jj) < 0 ) THEN
133                x_csthick(ji,jj) = ( z_lamda(ji,jj) * pp_v ) / z_fv(ji,jj)
134            ELSE
135                x_csthick(ji,jj) = 0.
136            ENDIF
137
138
139
140            ! Calculate the cool skin correction - only when the heat flux is out of the ocean
141            IF( tmask(ji,jj,1) == 1. .AND. x_csthick(ji,jj) /= 0. .AND. psqflux(ji,jj) < 0. ) THEN
142               x_csdsst(ji,jj) = ( psqflux(ji,jj) * x_csthick(ji,jj) ) / pp_k
143             ELSE
144               x_csdsst(ji,jj) = 0.
145            ENDIF
146
147         ENDDO
148      ENDDO
149
150   END SUBROUTINE diurnal_sst_coolskin_step
151
152
153END MODULE cool_skin
Note: See TracBrowser for help on using the repository browser.