Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
O
ocean-mixing
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Package registry
Model registry
Operate
Environments
Terraform modules
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
inmcm-mirror
ocean-mixing
Commits
ab23d4a8
Commit
ab23d4a8
authored
7 months ago
by
Evgeny Mortikov
Browse files
Options
Downloads
Patches
Plain Diff
state eq. implementation update
parent
1f3aca50
No related branches found
No related tags found
Loading
Changes
4
Show whitespace changes
Inline
Side-by-side
Showing
4 changed files
obl_bc.f90
+1
-1
1 addition, 1 deletion
obl_bc.f90
obl_main.f90
+1
-1
1 addition, 1 deletion
obl_main.f90
obl_scm.f90
+2
-2
2 additions, 2 deletions
obl_scm.f90
obl_state_eq.f90
+49
-30
49 additions, 30 deletions
obl_state_eq.f90
with
53 additions
and
34 deletions
obl_bc.f90
+
1
−
1
View file @
ab23d4a8
...
...
@@ -67,7 +67,7 @@ module obl_bc
theta_dyn
=
heat_flux
/
max
(
U_dyn
,
U_dyn_min
)
salin_dyn
=
salin_flux
/
max
(
U_dyn
,
U_dyn_min
)
rho_dyn
=
-
Rho_ref
*
(
alpha_state
*
theta_dyn
-
beta_state
*
salin_dyn
)
rho_dyn
=
Rho_ref
*
(
-
alpha_state
*
theta_dyn
+
beta_state
*
salin_dyn
)
end
subroutine
...
...
This diff is collapsed.
Click to expand it.
obl_main.f90
+
1
−
1
View file @
ab23d4a8
...
...
@@ -243,7 +243,7 @@ program obl_main
call
Theta_init
(
Theta
,
Theta_dev
,
grid
%
cz
,
grid
%
dz
)
!call Theta_init (Theta, nz, dz)
call
Salin_init
(
Salin
,
Salin_dev
,
grid
%
cz
,
grid
%
dz
)
call
solve_state_eq
(
Rho
,
Theta_dev
,
Salin_dev
,
grid
%
cz
)
call
get_density
(
Rho
,
Theta_dev
,
Salin_dev
,
grid
%
cz
)
call
U_init
(
U
,
grid
%
cz
)
call
V_init
(
V
,
grid
%
cz
)
...
...
This diff is collapsed.
Click to expand it.
obl_scm.f90
+
2
−
2
View file @
ab23d4a8
...
...
@@ -153,7 +153,7 @@ module obl_scm
real
,
intent
(
in
)
::
dt
! ----------------------------------------------------------------------------
call
solve_state_eq
(
Rho
,
Theta_dev
,
Salin_dev
,
grid
%
cz
)
call
get_density
(
Rho
,
Theta_dev
,
Salin_dev
,
grid
%
cz
)
end
subroutine
advance_state_eq_scm
...
...
This diff is collapsed.
Click to expand it.
obl_state_eq.f90
+
49
−
30
View file @
ab23d4a8
module
obl_state_eq
!< @brief state equation solver
! --------------------------------------------------------------------------------
! TO DO:
! -- differentiate between using 'dev' & 'full' values
! -- add non-linear state eq. support
! modules used
use
obl_math
! --------------------------------------------------------------------------------
! directives list
implicit
none
private
! public interface
! --------------------------------------------------------------------------------
public
::
get_density
,
get_density_dev
! --------------------------------------------------------------------------------
real
,
parameter
::
Rho_ref
=
999.1285
!< reference density, [kg / m**3]
!
real, parameter ::
Salin_0 = 35.0
!< reference
salin
ity, [
PSU
]
real
,
parameter
::
alpha_state
=
0.15063
/
Rho_ref
!< thermal expansion [1 /
K]
real
,
parameter
::
beta_state
=
7.8
*
0.0001
!< salinity expansion [1 /
PSU]
! --------------------------------------------------------------------------------
real
,
public
,
parameter
::
Rho_ref
=
999.1285
!< reference
dens
ity, [
kg / m**3
]
real
,
public
,
parameter
::
Theta_ref
=
288.15
!< reference temperature, [
K]
real
,
public
,
parameter
::
Salin_ref
=
35.0
!< reference salinity, [
PSU]
! reference temperature
real
,
parameter
::
Theta_ref
=
288.15
!< [K]
real
,
parameter
::
Salin_ref
=
35.0
!< reference salinity, [PSU]
!< linear state eq. parameters
real
,
public
,
parameter
::
alpha_state
=
0.15063
/
Rho_ref
!< thermal expansion [1 / K]
real
,
public
,
parameter
::
beta_state
=
7.8
*
0.0001
!< salinity expansion [1 / PSU]
! --------------------------------------------------------------------------------
contains
subroutine
solve_state_eq
(
Rho
,
Theta
,
Salin
,
nz
)
!< @brief state equation solver
real
,
intent
(
out
)
::
Rho
(
nz
)
!< density, [kg / m**3]
integer
,
intent
(
in
)
::
nz
!< number of z-steps
real
,
intent
(
in
)
::
Theta
(
nz
)
!< temperature, [K]
real
,
intent
(
in
)
::
Salin
(
nz
)
!< salinity, [PSU]
integer
::
k
!< counter
! --------------------------------------------------------------------------------
subroutine
get_density
(
Rho
,
Theta_dev
,
Salin_dev
,
cz
)
!< @brief state equation def.
! --------------------------------------------------------------------------------
integer
,
intent
(
in
)
::
cz
!< grid size
real
,
intent
(
out
)
::
Rho
(
cz
)
!< density, [kg / m**3]
real
,
intent
(
in
)
::
Theta_dev
(
cz
)
!< temperature, [K]
real
,
intent
(
in
)
::
Salin_dev
(
cz
)
!< salinity, [PSU]
do
k
=
1
,
nz
Rho
(
k
)
=
Rho_ref
*
(
1.0
-
alpha_state
*
(
Theta
(
k
))
+
beta_state
*
(
Salin
(
k
)))
integer
::
k
! --------------------------------------------------------------------------------
do
k
=
1
,
cz
Rho
(
k
)
=
Rho_ref
*
(
1.0
-
alpha_state
*
Theta_dev
(
k
)
+
beta_state
*
Salin_dev
(
k
))
end
do
end
subroutine
subroutine
solve_state_eq_dev
(
Rho_dev
,
Theta
,
Salin
,
nz
)
!< @brief state equation solver
! --------------------------------------------------------------------------------
subroutine
get_density_dev
(
Rho_dev
,
Theta_dev
,
Salin_dev
,
cz
)
!< @brief state equation def.
! --------------------------------------------------------------------------------
integer
,
intent
(
in
)
::
cz
!< grid size
real
,
intent
(
out
)
::
Rho_dev
(
nz
)
!< density, [kg / m**3]
integer
,
intent
(
in
)
::
nz
!< number of z-steps
real
,
intent
(
in
)
::
Theta
(
nz
)
!< temperature, [K]
real
,
intent
(
in
)
::
Salin
(
nz
)
!< salinity, [PSU]
integer
::
k
!< counter
real
,
intent
(
out
)
::
Rho_dev
(
cz
)
!< density, [kg / m**3]
real
,
intent
(
in
)
::
Theta_dev
(
cz
)
!< temperature, [K]
real
,
intent
(
in
)
::
Salin_dev
(
cz
)
!< salinity, [PSU]
integer
::
k
! --------------------------------------------------------------------------------
do
k
=
1
,
n
z
Rho_dev
(
k
)
=
Rho_ref
*
(
1.0
-
alpha_state
*
(
Theta
(
k
)
)
+
beta_state
*
(
Salin
(
k
))
)
-
Rho_ref
do
k
=
1
,
c
z
Rho_dev
(
k
)
=
Rho_ref
*
(
-
alpha_state
*
Theta
_dev
(
k
)
+
beta_state
*
Salin
_dev
(
k
))
end
do
end
subroutine
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment