Skip to content
Snippets Groups Projects
Commit f6d9f3ae authored by Evgeny Mortikov's avatar Evgeny Mortikov
Browse files

adding allocate/deallocate subroutines; I/O module

parent 747417f4
No related branches found
No related tags found
No related merge requests found
......@@ -2,6 +2,7 @@
gfortran -c -cpp -Wuninitialized srcF/sfx_phys_const.f90
gfortran -c -cpp -Wuninitialized srcF/sfx_common.f90
gfortran -c -cpp -Wuninitialized srcF/sfx_io.f90
gfortran -c -cpp -Wuninitialized srcF/sfx_data.f90
gfortran -c -cpp -Wuninitialized srcF/sfx_surface.f90
......@@ -13,5 +14,5 @@ gfortran -c -cpp -Wuninitialized srcF/sfx_esm_param.f90
gfortran -c -cpp -Wuninitialized srcF/sfx_esm.f90
gfortran -c -cpp -Wuninitialized srcF/sfx_main.f90
gfortran -o sfx.exe sfx_main.o sfx_log.o sfx_log_param.o sfx_esm.o sfx_esm_param.o sfx_surface.o sfx_data.o sfx_common.o sfx_phys_const.o
gfortran -o sfx.exe sfx_main.o sfx_log.o sfx_log_param.o sfx_esm.o sfx_esm_param.o sfx_surface.o sfx_data.o sfx_io.o sfx_common.o sfx_phys_const.o
rm *.o *.mod
......@@ -10,7 +10,7 @@ ifeq ($(COMPILER),gnu)
FC = gfortran
endif
OBJ_F90 = sfx_phys_const.o sfx_common.o sfx_data.o sfx_surface.o sfx_log_param.o sfx_log.o sfx_esm_param.o sfx_esm.o sfx_main.o
OBJ_F90 = sfx_phys_const.o sfx_common.o sfx_io.o sfx_data.o sfx_surface.o sfx_log_param.o sfx_log.o sfx_esm_param.o sfx_esm.o sfx_main.o
OBJ_F =
OBJ = $(OBJ_F90) $(OBJ_F)
......
module sfx_common
!> @brief surface flux code common subroutines
!> @brief surface flux model common subroutines
public
contains
......@@ -10,22 +10,28 @@ contains
! ----------------------------------------------------------------------------
implicit none
integer, intent(out) :: int
integer, intent(out) :: stat
integer, intent(out) :: stat !> output status, /= 0 signals ERROR
character(len = *), intent(in) :: str
! ----------------------------------------------------------------------------
read(str, * , iostat=stat) int
read(str, * , iostat = stat) int
end subroutine str2int
! ----------------------------------------------------------------------------
! ----------------------------------------------------------------------------
elemental function is_finite(value)
!> @brief check if value is finite
! ----------------------------------------------------------------------------
use ieee_arithmetic
implicit none
logical :: is_finite
real, intent(in) :: value
! ----------------------------------------------------------------------------
is_finite = ieee_is_finite(value)
end function is_finite
! ----------------------------------------------------------------------------
end module sfx_common
\ No newline at end of file
module sfx_data
!> @brief surface flux module data
!> @brief surface flux model module data
! modules used
! --------------------------------------------------------------------------------
......@@ -13,6 +13,9 @@ module sfx_data
! public interface
! --------------------------------------------------------------------------------
public :: allocate_meteo_vec, deallocate_meteo_vec
public :: allocate_sfx_vec, deallocate_sfx_vec
public :: push_sfx_data
! --------------------------------------------------------------------------------
......@@ -78,11 +81,94 @@ module sfx_data
contains
! --------------------------------------------------------------------------------
subroutine allocate_meteo_vec(meteo, n)
!> @brief allocate meteo data vector
! ----------------------------------------------------------------------------
type (meteoDataVecType), intent(inout) :: meteo
integer, intent(in) :: n
! ----------------------------------------------------------------------------
allocate(meteo%h(n))
allocate(meteo%U(n))
allocate(meteo%dT(n))
allocate(meteo%Tsemi(n))
allocate(meteo%dQ(n))
allocate(meteo%z0_m(n))
end subroutine allocate_meteo_vec
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
subroutine deallocate_meteo_vec(meteo)
!> @brief deallocate meteo data vector
! ----------------------------------------------------------------------------
type (meteoDataVecType), intent(inout) :: meteo
! ----------------------------------------------------------------------------
deallocate(meteo%h)
deallocate(meteo%U)
deallocate(meteo%dT)
deallocate(meteo%Tsemi)
deallocate(meteo%dQ)
deallocate(meteo%z0_m)
end subroutine deallocate_meteo_vec
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
subroutine allocate_sfx_vec(sfx, n)
!> @brief allocate surface fluxes data vector
! ----------------------------------------------------------------------------
type (sfxDataVecType), intent(inout) :: sfx
integer, intent(in) :: n
! ----------------------------------------------------------------------------
allocate(sfx%zeta(n))
allocate(sfx%Rib(n))
allocate(sfx%Re(n))
allocate(sfx%B(n))
allocate(sfx%z0_m(n))
allocate(sfx%z0_t(n))
allocate(sfx%Rib_conv_lim(n))
allocate(sfx%Cm(n))
allocate(sfx%Ct(n))
allocate(sfx%Km(n))
allocate(sfx%Pr_t_inv(n))
end subroutine allocate_sfx_vec
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
subroutine deallocate_sfx_vec(sfx)
!> @brief deallocate surface fluxes data vector
! ----------------------------------------------------------------------------
type (sfxDataVecType), intent(inout) :: sfx
! ----------------------------------------------------------------------------
deallocate(sfx%zeta)
deallocate(sfx%Rib)
deallocate(sfx%Re)
deallocate(sfx%B)
deallocate(sfx%z0_m)
deallocate(sfx%z0_t)
deallocate(sfx%Rib_conv_lim)
deallocate(sfx%Cm)
deallocate(sfx%Ct)
deallocate(sfx%Km)
deallocate(sfx%Pr_t_inv)
end subroutine deallocate_sfx_vec
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
subroutine push_sfx_data(sfx, sfx_cell, idx)
!> @brief helper subroutine for copying data in sfxDataVecType
! ----------------------------------------------------------------------------
type (sfxDataVecType), intent(inout) :: sfx
type (sfxDataType), intent(in) :: sfx_cell
integer, intent(in) :: idx
! ----------------------------------------------------------------------------
......
module sfx_io
!> @brief surface flux model I/O subroutines
implicit none
public
contains
!> @brief write data (2 vectors) in simple ascii format
! ----------------------------------------------------------------------------
subroutine write_ascii_vec2(fname, var1, var2, n, stat)
implicit none
integer, intent(out) :: stat
character(*), intent(in) :: fname
integer, intent(in) :: n
real, dimension(n), intent(in) :: var1, var2
! ----------------------------------------------------------------------------
! --- local variables
integer i
character(len = 7) str_stat
! ----------------------------------------------------------------------------
open(1, FILE = trim(fname), iostat = stat)
if (stat /= 0) return
do i = 1, n
write(1, *, iostat = stat) var1(i), var2(i)
if (stat /= 0) exit
end do
close(1, iostat = stat, STATUS = str_stat)
end subroutine write_ascii_vec2
! ----------------------------------------------------------------------------
!> @brief write data (11 vectors) in simple ascii format
! ----------------------------------------------------------------------------
subroutine write_ascii_vec11(fname, &
var1, var2, var3, var4, var5, var6, var7, var8, var9, var10, var11, &
n, fmt, stat)
implicit none
integer, intent(out) :: stat
character(*), intent(in) :: fname
character(*), intent(in) :: fmt
integer, intent(in) :: n
real, dimension(n), intent(in) :: var1, var2, var3, var4, var5
real, dimension(n), intent(in) :: var6, var7, var8, var9, var10, var11
! ----------------------------------------------------------------------------
! --- local variables
integer i
character(len = 7) str_stat
! ----------------------------------------------------------------------------
open(1, FILE = trim(fname), iostat = stat)
if (stat /= 0) return
do i = 1, n
write(1, fmt, iostat = stat) var1(i), var2(i), var3(i), var4(i), var5(i), &
var6(i), var7(i), var8(i), var9(i), var10(i), var11(i)
if (stat /= 0) exit
end do
close(1, iostat = stat, STATUS = str_stat)
end subroutine write_ascii_vec11
! ----------------------------------------------------------------------------
end module sfx_io
......@@ -4,6 +4,7 @@ program sfx_main
! --------------------------------------------------------------------------------
use sfx_phys_const
use sfx_common
use sfx_io
use sfx_data
use sfx_esm, only: &
......@@ -239,24 +240,8 @@ program sfx_main
!> @brief allocate input & output data
allocate(meteo%h(num))
allocate(meteo%U(num))
allocate(meteo%dT(num))
allocate(meteo%Tsemi(num))
allocate(meteo%dQ(num))
allocate(meteo%z0_m(num))
allocate(sfx%zeta(num))
allocate(sfx%Rib(num))
allocate(sfx%Re(num))
allocate(sfx%B(num))
allocate(sfx%z0_m(num))
allocate(sfx%z0_t(num))
allocate(sfx%Rib_conv_lim(num))
allocate(sfx%Cm(num))
allocate(sfx%Ct(num))
allocate(sfx%Km(num))
allocate(sfx%Pr_t_inv(num))
call allocate_meteo_vec(meteo, num)
call allocate_sfx_vec(sfx, num)
!> @brief read input data common parameters
......@@ -289,37 +274,19 @@ program sfx_main
!> @brief write output data
open(2, file = filename_out)
do i = 1, num
write(2, 20) sfx%zeta(i), sfx%Rib(i), &
sfx%Re(i), sfx%B(i), sfx%z0_m(i), sfx%z0_t(i), &
sfx%Rib_conv_lim(i), &
sfx%Cm(i),sfx%Ct(i), sfx%Km(i), sfx%Pr_t_inv(i)
enddo
close(2)
call write_ascii_vec11(filename_out, &
sfx%zeta, sfx%Rib, &
sfx%Re, sfx%B, sfx%z0_m, sfx%z0_t, &
sfx%Rib_conv_lim, &
sfx%Cm,sfx%Ct, sfx%Km, sfx%Pr_t_inv, num, '(11(f10.4,3x))', status)
!> @brief deallocate input & output data
deallocate(meteo%h)
deallocate(meteo%U)
deallocate(meteo%dT)
deallocate(meteo%Tsemi)
deallocate(meteo%dQ)
deallocate(meteo%z0_m)
deallocate(sfx%zeta)
deallocate(sfx%Rib)
deallocate(sfx%Re)
deallocate(sfx%B)
deallocate(sfx%z0_m)
deallocate(sfx%z0_t)
deallocate(sfx%Rib_conv_lim)
deallocate(sfx%Cm)
deallocate(sfx%Ct)
deallocate(sfx%Km)
deallocate(sfx%Pr_t_inv)
! *: remove format(10) if not needed
call deallocate_meteo_vec(meteo)
call deallocate_sfx_vec(sfx)
! *: remove formats: not needed
10 format (f8.4,2x,f8.4)
20 format (11(f10.4,3x))
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment