Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
S
sfx
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
sfx
Commits
0f9143aa
Commit
0f9143aa
authored
2 weeks ago
by
Debolskiy Andrey
Browse files
Options
Downloads
Patches
Plain Diff
introduce lake model api (needed since it uses double precision)
parent
1b0d54db
No related branches found
No related tags found
No related merge requests found
Pipeline
#1803
canceled
2 weeks ago
Stage: build
Stage: test
Changes
3
Pipelines
1
Hide whitespace changes
Inline
Side-by-side
Showing
3 changed files
CMakeLists.txt
+1
-0
1 addition, 0 deletions
CMakeLists.txt
srcF/sfx_api_lake.f90
+77
-0
77 additions, 0 deletions
srcF/sfx_api_lake.f90
srcF/sfx_surface.f90
+2
-0
2 additions, 0 deletions
srcF/sfx_surface.f90
with
80 additions
and
0 deletions
CMakeLists.txt
+
1
−
0
View file @
0f9143aa
...
...
@@ -100,6 +100,7 @@ set(SOURCES_F
srcF/sfx_fc_wrapper.F90
srcF/sfx_api_inmcm.f90
srcF/sfx_api_term.f90
srcF/sfx_api_lake.f90
)
set
(
DIAG
...
...
This diff is collapsed.
Click to expand it.
srcF/sfx_api_lake.f90
0 → 100644
+
77
−
0
View file @
0f9143aa
!> @brief sfx-inmcm coupling API
module
sfx_api_lake
! modules used
! --------------------------------------------------------------------------------
use
sfx_data
! --------------------------------------------------------------------------------
! directives list
! --------------------------------------------------------------------------------
implicit
none
private
! --------------------------------------------------------------------------------
! public interface
! --------------------------------------------------------------------------------
public
::
lake_to_sfx_in_cell
,
sfx_to_lake_out_cell
! --------------------------------------------------------------------------------
contains
! --------------------------------------------------------------------------------
subroutine
lake_to_sfx_in_cell
(
meteo
,
arg
,
IVEG_sfx
,
depth_inm
,
lai_inm
)
!> @brief converts legacy arg [AR1 INMCM format but double pres] array to sfx meteo input
! ----------------------------------------------------------------------------
use
,
intrinsic
::
iso_c_binding
,
only
:
real_d
=>
c_double
,
&
! 8-byte real
real_f
=>
c_float
! 4-byte real
implicit
none
type
(
meteoDataType
),
intent
(
inout
)
::
meteo
real
(
kind
=
real_d
),
dimension
(
6
),
intent
(
in
)
::
arg
integer
,
intent
(
in
)
::
IVEG_sfx
real
(
kind
=
real_d
),
intent
(
in
)
::
depth_inm
real
(
kind
=
real_d
),
intent
(
in
)
::
lai_inm
! ----------------------------------------------------------------------------
meteo
%
U
=
real
(
arg
(
1
),
real_f
)
meteo
%
dT
=
real
(
arg
(
2
),
real_f
)
meteo
%
Tsemi
=
real
(
arg
(
3
),
real_f
)
meteo
%
dQ
=
real
(
arg
(
4
),
real_f
)
meteo
%
h
=
real
(
arg
(
5
),
real_f
)
meteo
%
z0_m
=
real
(
arg
(
6
),
real_f
)
meteo
%
depth
=
real
(
depth_inm
,
real_f
)
meteo
%
lai
=
real
(
lai_inm
,
real_f
)
meteo
%
surface_type
=
IVEG_sfx
!write(*,*) 'surface_type, IVEG_sfx', meteo%surface_type, IVEG_sfx
end
subroutine
lake_to_sfx_in_cell
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
subroutine
sfx_to_lake_out_cell
(
arg
,
sfx
)
!> @brief converts sfx cell output to legacy arg [AR2 INMCM format] array
! ----------------------------------------------------------------------------
use
,
intrinsic
::
iso_c_binding
,
only
:
real_d
=>
c_double
! 8-byte real
implicit
none
type
(
sfxDataType
),
intent
(
in
)
::
sfx
real
(
kind
=
real_d
),
dimension
(
11
),
intent
(
inout
)
::
arg
! ----------------------------------------------------------------------------
arg
(
1
)
=
dble
(
sfx
%
zeta
)
arg
(
2
)
=
dble
(
sfx
%
Rib
)
arg
(
3
)
=
dble
(
sfx
%
Re
)
arg
(
4
)
=
dble
(
sfx
%
B
)
arg
(
5
)
=
dble
(
sfx
%
z0_m
)
arg
(
6
)
=
dble
(
sfx
%
z0_t
)
!arg(7) = 0.0 ! arg(7) is never used in legacy code
arg
(
8
)
=
dble
(
sfx
%
Cm
)
arg
(
9
)
=
dble
(
sfx
%
Ct
)
arg
(
10
)
=
dble
(
sfx
%
Km
)
arg
(
11
)
=
dble
(
sfx
%
Pr_t_inv
)
end
subroutine
sfx_to_lake_out_cell
! --------------------------------------------------------------------------------
end
module
sfx_api_lake
This diff is collapsed.
Click to expand it.
srcF/sfx_surface.f90
+
2
−
0
View file @
0f9143aa
...
...
@@ -322,6 +322,8 @@ contains
#if defined(INCLUDE_CXX)
subroutine
set_c_struct_sfx_surface_param_values
(
surface_param
)
use
sfx_data
use
sfx_z0m_all_surface
use
sfx_z0t_all_surface
implicit
none
type
(
sfx_surface_param
),
intent
(
inout
)
::
surface_param
surface_param
%
surface_ocean
=
surface_ocean
...
...
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