Line data Source code
1 : module co2_data_flux
2 :
3 : !-------------------------------------------------------------------------------
4 : ! utilities for reading and interpolating co2 surface fluxes
5 : !-------------------------------------------------------------------------------
6 :
7 : use shr_kind_mod, only: r8 => shr_kind_r8, cl => shr_kind_cl
8 : use input_data_utils, only: time_coordinate
9 : use cam_abortutils, only: endrun
10 :
11 : implicit none
12 : private
13 : save
14 :
15 : ! Public interfaces
16 : public co2_data_flux_type
17 :
18 : public co2_data_flux_init
19 : public co2_data_flux_advance
20 :
21 : !-------------------------------------------------------------------------------
22 :
23 : type :: co2_data_flux_type
24 : character(len=cl) :: filename
25 : character(len=cl) :: varname
26 : logical :: initialized
27 : type(time_coordinate) :: time_coord
28 : real(r8), pointer :: co2bdy(:,:,:) ! bracketing data (pcols,begchunk:endchunk,2)
29 : real(r8), pointer :: co2flx(:,:) ! Interpolated output (pcols,begchunk:endchunk)
30 : end type co2_data_flux_type
31 :
32 : ! dimension names for physics grid (physgrid)
33 : logical :: dimnames_set = .false.
34 : character(len=8) :: dim1name, dim2name
35 :
36 : !===============================================================================
37 : contains
38 : !===============================================================================
39 :
40 0 : subroutine co2_data_flux_init (input_file, varname, xin)
41 :
42 : !-------------------------------------------------------------------------------
43 : ! Initialize co2_data_flux_type instance
44 : ! including initial read of input and interpolation to the current timestep
45 : !-------------------------------------------------------------------------------
46 :
47 : use ioFileMod, only: getfil
48 : use ppgrid, only: begchunk, endchunk, pcols
49 : use cam_grid_support, only: cam_grid_id, cam_grid_check
50 : use cam_grid_support, only: cam_grid_get_dim_names
51 :
52 : ! Arguments
53 : character(len=*), intent(in) :: input_file
54 : character(len=*), intent(in) :: varname
55 : type(co2_data_flux_type), intent(inout) :: xin
56 :
57 : ! Local variables
58 : integer :: grid_id
59 : real(r8) :: dtime
60 : character(len=*), parameter :: subname = 'co2_data_flux_init'
61 : !----------------------------------------------------------------------------
62 :
63 0 : if (.not. dimnames_set) then
64 0 : grid_id = cam_grid_id('physgrid')
65 0 : if (.not. cam_grid_check(grid_id)) then
66 0 : call endrun(subname // ': ERROR: no "physgrid" grid')
67 : endif
68 0 : call cam_grid_get_dim_names(grid_id, dim1name, dim2name)
69 0 : dimnames_set = .true.
70 : end if
71 :
72 0 : call getfil(input_file, xin%filename)
73 0 : xin%varname = varname
74 0 : xin%initialized = .false.
75 :
76 0 : dtime = 1.0_r8 - 200.0_r8 / 86400.0_r8
77 0 : call xin%time_coord%initialize(input_file, delta_days=dtime)
78 :
79 : allocate( xin%co2bdy(pcols,begchunk:endchunk,2), &
80 0 : xin%co2flx(pcols,begchunk:endchunk) )
81 :
82 0 : call co2_data_flux_advance(xin)
83 :
84 0 : xin%initialized = .true.
85 :
86 0 : end subroutine co2_data_flux_init
87 :
88 : !===============================================================================
89 :
90 0 : subroutine co2_data_flux_advance (xin)
91 :
92 : !-------------------------------------------------------------------------------
93 : ! Advance the contents of a co2_data_flux_type instance
94 : ! including reading new data, if necessary
95 : !-------------------------------------------------------------------------------
96 :
97 0 : use cam_pio_utils, only: cam_pio_openfile
98 : use ncdio_atm, only: infld
99 : use pio, only: file_desc_t, pio_nowrite, pio_closefile
100 : use ppgrid, only: begchunk, endchunk, pcols
101 :
102 : ! Arguments
103 : type(co2_data_flux_type), intent(inout) :: xin
104 :
105 : ! Local variables
106 : character(len=*), parameter :: subname = 'co2_data_flux_advance'
107 : logical :: read_data
108 : integer :: indx2_pre_adv
109 : type(file_desc_t) :: fh_co2_data_flux
110 : logical :: found
111 :
112 : !----------------------------------------------------------------------------
113 :
114 0 : read_data = xin%time_coord%read_more() .or. .not. xin%initialized
115 :
116 0 : indx2_pre_adv = xin%time_coord%indxs(2)
117 :
118 0 : call xin%time_coord%advance()
119 :
120 0 : if ( read_data ) then
121 :
122 0 : call cam_pio_openfile(fh_co2_data_flux, trim(xin%filename), PIO_NOWRITE)
123 :
124 : ! read time-level 1
125 : ! skip the read if the needed vals are present in time-level 2
126 0 : if (xin%initialized .and. xin%time_coord%indxs(1) == indx2_pre_adv) then
127 0 : xin%co2bdy(:,:,1) = xin%co2bdy(:,:,2)
128 : else
129 : call infld(trim(xin%varname), fh_co2_data_flux, dim1name, dim2name, &
130 0 : 1, pcols, begchunk, endchunk, xin%co2bdy(:,:,1), found, &
131 0 : gridname='physgrid', timelevel=xin%time_coord%indxs(1))
132 0 : if (.not. found) then
133 0 : call endrun(subname // ': ERROR: ' // trim(xin%varname) // ' not found')
134 : endif
135 : endif
136 :
137 : ! read time-level 2
138 : call infld(trim(xin%varname), fh_co2_data_flux, dim1name, dim2name, &
139 0 : 1, pcols, begchunk, endchunk, xin%co2bdy(:,:,2), found, &
140 0 : gridname='physgrid', timelevel=xin%time_coord%indxs(2))
141 0 : if (.not. found) then
142 0 : call endrun(subname // ': ERROR: ' // trim(xin%varname) // ' not found')
143 : endif
144 :
145 0 : call pio_closefile(fh_co2_data_flux)
146 : endif
147 :
148 : ! interpolate between time-levels
149 : ! If time:bounds is in the dataset, and the dataset calendar is compatible with CAM's,
150 : ! then the time_coordinate class will produce time_coord%wghts(2) == 0.0,
151 : ! generating fluxes that are piecewise constant in time.
152 :
153 0 : if (xin%time_coord%wghts(2) == 0.0_r8) then
154 0 : xin%co2flx(:,:) = xin%co2bdy(:,:,1)
155 : else
156 0 : xin%co2flx(:,:) = xin%co2bdy(:,:,1) + &
157 0 : xin%time_coord%wghts(2) * (xin%co2bdy(:,:,2) - xin%co2bdy(:,:,1))
158 : endif
159 :
160 0 : end subroutine co2_data_flux_advance
161 :
162 : !===============================================================================
163 :
164 0 : end module co2_data_flux
|