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/branches/2019/dev_r11943_MERGE_2019/src/OCE/DIU – NEMO

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/DIU/diu_coolskin.F90 @ 11949

Last change on this file since 11949 was 11949, checked in by acc, 4 years ago

Merge in changes from 2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps. This just creates a fresh copy of this branch to use as the merge base. See ticket #2341

  • Property svn:keywords set to Id
File size: 6.6 KB
RevLine 
[10989]1MODULE diu_coolskin
[5676]2   !!======================================================================
[10989]3   !!                    ***  MODULE  diu_coolskin  ***
[5676]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   !!----------------------------------------------------------------------
[7646]12   !!   diurnal_sst_coolskin_init : initialisation of the cool skin
13   !!   diurnal_sst_coolskin_step : time-stepping  of the cool skin corrections
[5676]14   !!----------------------------------------------------------------------
15   USE par_kind
16   USE phycst
17   USE dom_oce
18   USE in_out_manager
19   USE sbc_oce
[6493]20   USE lib_mpp
[5676]21   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
22   
23   IMPLICIT NONE
[7646]24   PRIVATE
25
[5676]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 "vectopt_loop_substitute.h90"
[7646]44   !!----------------------------------------------------------------------
[10068]45   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
46   !! $Id$
47   !! Software governed by the CeCILL license (see ./LICENSE)
[7646]48   !!----------------------------------------------------------------------   
[5676]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.
[7646]65      !
[5676]66   END SUBROUTINE diurnal_sst_coolskin_init
[7646]67
68
[5676]69   SUBROUTINE diurnal_sst_coolskin_step(psqflux, pstauflux, psrho, rdt)
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) :: 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
[7646]95      !!----------------------------------------------------------------------
96      !
[10989]97      IF( .NOT. ln_blk )   CALL ctl_stop("diu_coolskin.f90: diurnal flux processing only implemented for bulk forcing")
[7646]98      !
[5676]99      DO jj = 1,jpj
100         DO ji = 1,jpi
[7646]101            !
[5676]102            ! Calcualte wind speed from wind stress and friction velocity
103            IF( tmask(ji,jj,1) == 1. .AND. pstauflux(ji,jj) /= 0 .AND. psrho(ji,jj) /=0 ) THEN
104               z_fv(ji,jj) = SQRT( pstauflux(ji,jj) / psrho(ji,jj) )
105               z_wspd(ji,jj) = SQRT( pstauflux(ji,jj) / ( pp_cda * pp_rhoa ) )
106            ELSE
107               z_fv(ji,jj) = 0.
108               z_wspd(ji,jj) = 0.     
109            ENDIF
[7646]110            !
[5676]111            ! Calculate gamma function which is dependent upon wind speed
112            IF( tmask(ji,jj,1) == 1. ) THEN
113               IF( ( z_wspd(ji,jj) <= 7.5 ) ) z_gamma(ji,jj) = ( 0.2 * z_wspd(ji,jj) ) + 0.5
114               IF( ( z_wspd(ji,jj) > 7.5 ) .AND. ( z_wspd(ji,jj) < 10. ) ) z_gamma(ji,jj) = ( 1.6 * z_wspd(ji,jj) ) - 10.
115               IF( ( z_wspd(ji,jj) >= 10. ) ) z_gamma(ji,jj) = 6.
116            ENDIF
[7646]117            !
[5676]118            ! Calculate 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
[7646]124            !
[5676]125            ! Calculate the cool skin thickness - only when heat flux is out of the ocean
126            IF( tmask(ji,jj,1) == 1. .AND. z_fv(ji,jj) /= 0 .AND. psqflux(ji,jj) < 0 ) THEN
[7646]127               x_csthick(ji,jj) = ( z_lamda(ji,jj) * pp_v ) / z_fv(ji,jj)
[5676]128            ELSE
[7646]129               x_csthick(ji,jj) = 0.
[5676]130            ENDIF
[7646]131            !
[5676]132            ! Calculate the cool skin correction - only when the heat flux is out of the ocean
133            IF( tmask(ji,jj,1) == 1. .AND. x_csthick(ji,jj) /= 0. .AND. psqflux(ji,jj) < 0. ) THEN
134               x_csdsst(ji,jj) = ( psqflux(ji,jj) * x_csthick(ji,jj) ) / pp_k
135             ELSE
136               x_csdsst(ji,jj) = 0.
137            ENDIF
[7646]138            !
139         END DO
140      END DO
141      !
[5676]142   END SUBROUTINE diurnal_sst_coolskin_step
143
[7646]144   !!======================================================================
[10989]145END MODULE diu_coolskin
Note: See TracBrowser for help on using the repository browser.