-
Notifications
You must be signed in to change notification settings - Fork 119
/
Copy pathcubed_sphere_inc_mod.F90
66 lines (52 loc) · 2.52 KB
/
cubed_sphere_inc_mod.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
module cubed_sphere_inc_mod
use tracer_manager_mod, only: get_tracer_index, get_number_tracers, get_tracer_names
use field_manager_mod, only: MODEL_ATMOS
use fv_arrays_mod, only: fv_atmos_type
use fms2_io_mod, only: open_file, close_file, read_data, variable_exists, &
FmsNetcdfDomainFile_t, register_axis
implicit none
type increment_data_type
real, allocatable :: ua_inc(:,:,:)
real, allocatable :: va_inc(:,:,:)
real, allocatable :: temp_inc(:,:,:)
real, allocatable :: delp_inc(:,:,:)
real, allocatable :: delz_inc(:,:,:)
real, allocatable :: tracer_inc(:,:,:,:)
end type increment_data_type
public read_cubed_sphere_inc, increment_data_type
contains
!----------------------------------------------------------------------------------------
subroutine read_cubed_sphere_inc(fname, increment_data, Atm)
character(*), intent(in) :: fname
type(increment_data_type), intent(inout) :: increment_data
type(fv_atmos_type), intent(in) :: Atm
type(FmsNetcdfDomainFile_t) :: fileobj
integer :: itracer, ntracers
character(len=64) :: tracer_name
! Get various dimensions
call get_number_tracers(MODEL_ATMOS, num_tracers=ntracers)
! Open file
if ( open_file(fileobj, trim(fname), 'read', Atm%domain) ) then
! Register axes
call register_axis(fileobj, 'xaxis_1', 'x')
call register_axis(fileobj, 'yaxis_1', 'y')
! Read increments
call read_data(fileobj, 'u_inc', increment_data%ua_inc)
call read_data(fileobj, 'v_inc', increment_data%va_inc)
call read_data(fileobj, 'T_inc', increment_data%temp_inc)
call read_data(fileobj, 'delp_inc', increment_data%delp_inc)
if ( .not. Atm%flagstruct%hydrostatic ) then
call read_data(fileobj, 'delz_inc', increment_data%delz_inc)
end if
! Read tracer increments
do itracer = 1,ntracers
call get_tracer_names(MODEL_ATMOS, itracer, tracer_name)
if ( variable_exists(fileobj, trim(tracer_name)//'_inc') ) then
call read_data(fileobj, trim(tracer_name)//'_inc', increment_data%tracer_inc(:,:,:,itracer))
end if
end do
call close_file(fileobj)
end if
end subroutine read_cubed_sphere_inc
!----------------------------------------------------------------------------------------
end module cubed_sphere_inc_mod