SOCA
soca_kst_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2017-2021 UCAR
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 
6 !> variable transform: S/T balance
8 
9 use kinds, only: kind_real
10 use soca_utils, only: soca_diff
11 
12 implicit none
13 private
14 public :: soca_soft_jacobian
15 
16 
17 !> Hold the configuration and jacobians for the s/t balance transform
18 !!
19 !! should be populated by calls to soca_kst_mod::soca_kst::soca_soft_jacobian
20 !! \see soca_balance_mod::soca_balance_setup
21 type, public :: soca_kst
22  real(kind=kind_real) :: dsdtmax !< 1.0 [psu/K]
23  real(kind=kind_real) :: dsdzmin !< 3.0e-3 [psu/m]
24  real(kind=kind_real) :: dtdzmin !< 1.0e-3 [K/m]
25  integer :: nlayers
26  real(kind=kind_real), allocatable :: jacobian(:,:,:) !< dS/dT(i,j,k)
27 end type soca_kst
28 
29 
30 ! ------------------------------------------------------------------------------
31 contains
32 ! ------------------------------------------------------------------------------
33 
34 
35 ! ------------------------------------------------------------------------------
36 !> Jacobian of Sb=S(T) relative to the reference state t, s. jac=dS/dT at (t,s)
37 !!
38 !! \relates soca_kst_mod::soca_kst
39 !! \param s: Background practical salinity [g/kg]
40 !! \param t: Background potential Temperature [deg C]
41 !! \param h: Layer thickness [m]
42 !! \param jac: Jacobian [ds1/dt1, ...,dsN/dtN]; [psu/deg C]
43 subroutine soca_soft_jacobian (jac, t, s, h, dsdtmax, dsdzmin, dtdzmin)
44  real(kind=kind_real), intent(in) :: t(:), s(:), h(:)
45  real(kind=kind_real), intent(in) :: dsdtmax, dsdzmin, dtdzmin
46  real(kind=kind_real), allocatable, intent(inout) :: jac(:) ! jac=ds/dt
47 
48  real(kind=kind_real), allocatable :: dtdz(:), dsdz(:)
49  integer :: nl, z
50  real(kind=kind_real) :: j
51 
52  ! Allocate
53  nl = size(t,1)
54  allocate(dtdz(nl),dsdz(nl))
55  if (.not.allocated(jac)) allocate(jac(nl))
56 
57  call soca_diff(dtdz,t,h)
58  call soca_diff(dsdz,s,h)
59 
60  jac = 0.0
61  do z=1,nl
62  jac(z) = 0.0
63 
64  ! Limit application of soft according to configuration
65  if ( abs(dtdz(z)) < dtdzmin ) cycle
66  if ( abs(dsdz(z)) < dsdzmin ) cycle
67 
68  ! Jacobian of soft
69  j=dsdz(z)/dtdz(z)
70 
71  ! Limit application of soft according to configuration
72  if ( abs(j) > dsdtmax ) cycle
73 
74  ! if we reach this point in the code, the jacobian is usable
75  jac(z) = j;
76  end do
77 
78 end subroutine soca_soft_jacobian
79 
80 end module soca_kst_mod
variable transform: S/T balance
Definition: soca_kst_mod.F90:7
various utility functions
Definition: soca_utils.F90:7
subroutine, public soca_diff(dvdz, v, h)
Definition: soca_utils.F90:96
Hold the configuration and jacobians for the s/t balance transform.