REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_NCTimeSeriesRiver.cpp
Go to the documentation of this file.
2#include "REMORA_NCFile.H"
3
4#include "AMReX_ParallelDescriptor.H"
5
6#include <string>
7
8#ifdef REMORA_USE_NETCDF
9
10/**
11 * @param[in ] a_file_names vector of file name(s) to read from
12 * @param[in ] a_field_name name of field to read in
13 * @param[in ] a_time_name name of time variable in NetCDF file
14 * @param[in ] a_nz number of vertical levels in domain
15 * @param[in ] a_use_vert_integ whether the data in the file is vertically integrated
16 */
17NCTimeSeriesRiver::NCTimeSeriesRiver (const amrex::Vector<std::string>& a_file_names, const std::string a_field_name,
18 const std::string a_time_name,
19 const int a_nz, const int a_use_vert_integ) {
20 file_names.assign(a_file_names.begin(), a_file_names.end());
21 time_name = a_time_name;
22 field_name = a_field_name;
23 nz = a_nz;
24 use_vert_integ = a_use_vert_integ;
25}
26
28 // open file
29 amrex::Print() << "Loading " << field_name << " from rivers NetCDF file(s)" << std::endl;
30
31 // The time field can have any number of names, depending on the field.
32 // If not specified in input file (time_name.empty()) then set it by default
33 if (time_name.empty())
34 {
35 time_name = "river_time";
36 }
37
38 for (int ifile = 0; ifile < file_names.size(); ++ifile) {
39 const std::string& file_name = file_names[ifile];
40
41 // Check units of time stamps; should be days
42 std::string unit_str = ReadNetCDFVarAttrStr(file_name, time_name, "units"); // works on proc 0
43 if (amrex::ParallelDescriptor::IOProcessor())
44 {
45 if (unit_str.find("days") == std::string::npos) {
46 amrex::Print() << "Units of river_time given as: " << unit_str << std::endl;
47 amrex::Abort("Units must be in days.");
48 }
49 }
50
51 // get times and put in array
52 using RARRAY = NDArray<amrex::Real>;
53 amrex::Vector<RARRAY> array_ts(1);
54 ReadNetCDFFile(file_name, {time_name}, array_ts); // filled only on proc 0
55 if (amrex::ParallelDescriptor::IOProcessor())
56 {
57 int ntimes_io = array_ts[0].get_vshape()[0];
58 for (int nt(0); nt < ntimes_io; nt++)
59 {
60 // Convert river time from days to seconds
61 river_times.push_back((*(array_ts[0].get_data() + nt)) * amrex::Real(60.0) * amrex::Real(60.0) * amrex::Real(24.0));
62 file_for_time.push_back(ifile);
63 file_itime_offset.push_back(nt);
64 }
65 }
66 }
67 int ntimes = river_times.size();
68 // Only do checks on IO processors since river_times isn't populated on other ranks yet
69 if (amrex::ParallelDescriptor::IOProcessor()) {
70 AMREX_ASSERT(std::is_sorted(river_times.begin(), river_times.end()));
71 if (ntimes <= 1) {
72 amrex::Error("River data must be given at at least two times");
73 }
74 }
75 int ioproc = amrex::ParallelDescriptor::IOProcessorNumber();
76 amrex::ParallelDescriptor::Bcast(&ntimes,1,ioproc);
77 if (!(amrex::ParallelDescriptor::IOProcessor())) {
78 river_times.resize(ntimes);
79 file_for_time.resize(ntimes);
80 file_itime_offset.resize(ntimes);
81 }
82 amrex::ParallelDescriptor::Bcast(river_times.data(), river_times.size(), ioproc);
83 amrex::ParallelDescriptor::Bcast(file_for_time.data(), file_for_time.size(), ioproc);
84 amrex::ParallelDescriptor::Bcast(file_itime_offset.data(), file_itime_offset.size(), ioproc);
85
86 auto ncf = ncutils::NCFile::open(file_names[0], NC_NOCLOBBER );
87 if (amrex::ParallelDescriptor::IOProcessor())
88 {
89 std::vector<MPI_Offset> shape = ncf.var(field_name).shape();
90 if (shape.size() == 2) {
91 has_z = 0;
92 nriv = shape[1];
93 } else if (shape.size() == 3) {
94 has_z = 1;
95 nriv = shape[2];
96 } else {
97 amrex::Abort("River field shape not 2 or 3");
98 }
99 }
100 amrex::ParallelDescriptor::Bcast(&has_z, 1, ioproc);
101 amrex::ParallelDescriptor::Bcast(&nriv, 1, ioproc);
102
103 amrex::Box vshape_box(amrex::IntVect(0,0,0),amrex::IntVect(nriv,0,nz));
104 if (!has_z && !use_vert_integ) {
105 amrex::Vector<amrex::FArrayBox*> NC_fabs;
106 amrex::Vector<std::string> NC_names;
107 amrex::Vector<enum NC_Data_Dims_Type> NC_dim_types;
108 amrex::Print() << "Reading in river_Vshape from " << file_names[0] << std::endl;
109
110 fab_vshape = new amrex::FArrayBox();
111 NC_fabs.push_back(fab_vshape); NC_names.push_back("river_Vshape");
112 NC_dim_types.push_back(NC_Data_Dims_Type::BT_Riv);
113 BuildFABsFromNetCDFFile<amrex::FArrayBox,amrex::Real>(vshape_box, file_names[0], NC_names, NC_dim_types, NC_fabs);
114 }
115
116 nzbox = (use_vert_integ) ? 1 : nz;
117 amrex::Box riv_box(amrex::IntVect(0,0,0),amrex::IntVect(nriv,0,nzbox));
118#ifdef AMREX_USE_GPU
119 // It's possible there should be a different arena
120 fab_before = new amrex::FArrayBox(riv_box,1,amrex::The_Pinned_Arena());
121 fab_after = new amrex::FArrayBox(riv_box,1,amrex::The_Pinned_Arena());
122 fab_interp = new amrex::FArrayBox(riv_box,1,amrex::The_Pinned_Arena());
123#else
124 fab_before = new amrex::FArrayBox(riv_box,1);
125 fab_after = new amrex::FArrayBox(riv_box,1);
126 fab_interp = new amrex::FArrayBox(riv_box,1);
127#endif
128
129 // dummy initialization
130 i_time_before = -100;
131}
132
134 // Figure out time index:
135 AMREX_ASSERT(time >= river_times[0]);
136 AMREX_ASSERT(time <= river_times[river_times.size()-1]);
137 int i_time_before_old = i_time_before;
138 for (int nt=0; nt < river_times.size()-1; nt++) {
139 if ((river_times[nt] <= time) and (river_times[nt+1] >= time)) {
140 i_time_before = nt;
142 time_after = river_times[nt+1];
143 break;
144 }
145 }
146
147 int i_time_after = i_time_before + 1;
148 if (i_time_before_old + 1 == i_time_before) {
149 // swap data vectors so we only have to read in one MultiFab
150 std::swap(fab_before, fab_after);
151 read_in_at_time(fab_after, i_time_after);
152 } else if (i_time_before_old != i_time_before) {
153 read_in_at_time(fab_after, i_time_after);
155 }
156
157 amrex::Real dt = time_after - time_before;
158 amrex::Real time_before_copy = time_before;
159
160 amrex::Box fab_domain(amrex::IntVect(0,0,0), amrex::IntVect(nriv-1,0,nzbox-1));
161 auto interp_array = fab_interp->array();
162 auto before_array = fab_before->array();
163 auto after_array = fab_after->array();
164 amrex::ParallelFor(fab_domain, [=] AMREX_GPU_DEVICE (int r, int , int k) {
165 interp_array(r,0,k) = before_array(r,0,k) + (time - time_before_copy) * (after_array(r,0,k) - before_array(r,0,k)) / dt;
166 });
167}
168
169void NCTimeSeriesRiver::read_in_at_time (amrex::FArrayBox* fab_dat, int itime) {
170 amrex::FArrayBox NC_fab;
171 amrex::Vector<amrex::FArrayBox*> NC_fabs;
172 amrex::Vector<std::string> NC_names;
173 amrex::Vector<enum NC_Data_Dims_Type> NC_dim_types;
174 // actual dims don't really matter here; only lower is used in call
175
176 const std::string& file_name = file_names[file_for_time[itime]];
177 const int itime_offset = file_itime_offset[itime];
178
179 amrex::Print() << "Reading in " << field_name << " at time index " << itime
180 << " from " << file_name << std::endl;
181
182 NC_fabs.push_back(&NC_fab) ; NC_names.push_back(field_name);
183
184 if (has_z) {
185 NC_dim_types.push_back(NC_Data_Dims_Type::Time_BT_Riv);
186 } else {
187 NC_dim_types.push_back(NC_Data_Dims_Type::Time_Riv);
188 }
189
190 amrex::Box riv_domain(amrex::IntVect(0,0,0), amrex::IntVect(nriv-1,0,nz-1));
191 amrex::Box fab_domain(amrex::IntVect(0,0,0), amrex::IntVect(nriv-1,0,nzbox-1));
192
193 BuildFABsFromNetCDFFile<amrex::FArrayBox,amrex::Real>(riv_domain, file_name, NC_names, NC_dim_types, NC_fabs, true, itime_offset);
194
195 auto dat_array = fab_dat->array();
196 auto tmp_array = NC_fabs[0]->const_array();
197 if (has_z || use_vert_integ) {
198 amrex::ParallelFor(fab_domain, [=] AMREX_GPU_DEVICE (int r, int , int k) {
199 dat_array(r,0,k) = tmp_array(r,0,k);
200 });
201 } else {
202 auto array_vshape = fab_vshape->const_array();
203 amrex::ParallelFor(fab_domain, [=] AMREX_GPU_DEVICE (int r, int , int k) {
204 dat_array(r,0,k) = tmp_array(r,0,0) * array_vshape(r,0,k);
205 });
206 }
207}
208#endif // REMORA_USE_NETCDF
void ReadNetCDFFile(const std::string &fname, amrex::Vector< std::string > names, amrex::Vector< NDArray< DType > > &arrays, bool one_time=false, int fill_time=0)
Read in data from netcdf file and save to data arrays.
std::string ReadNetCDFVarAttrStr(const std::string &fname, const std::string &var_name, const std::string &attr_name)
Helper function for reading a single variable attribute.
int nriv
Number of rivers.
amrex::Real time_before
Time in ocean_time immediately before the last time interpolated to.
int nzbox
How many vertical cells there are in the data fabs.
void update_interpolated_to_time(amrex::Real time)
Calculate.
amrex::Vector< int > file_itime_offset
Offset to access a particular time within its file.
amrex::Vector< std::string > file_names
File names to read from.
std::string time_name
Field name for time series in netcdf file.
amrex::Vector< amrex::Real > river_times
Time points in netcdf file.
NCTimeSeriesRiver(const amrex::Vector< std::string > &a_file_names, const std::string a_field_name, const std::string a_time_name, const int a_nz, const int a_use_vert_integ=0)
std::string field_name
Field name in netcdf file.
amrex::FArrayBox * fab_vshape
Vshape data if needed.
amrex::FArrayBox * fab_after
int i_time_before
Time index immediately before the last time interpolated to.
void read_in_at_time(amrex::FArrayBox *vec, int itime)
amrex::FArrayBox * fab_before
FABs to pointers of river data.
int nz
Number of vertical points.
amrex::Vector< int > file_for_time
File index to access a particular time.
amrex::FArrayBox * fab_interp
Container for interpolated data; Only used if save_interpolated == true.
int has_z
Whether the field is specified in the z-dimension.
amrex::Real time_after
Time in ocean_time immediately after the last time interpolated to.
static NCFile open(const std::string &name, const int cmode=NC_NOWRITE, MPI_Comm comm=MPI_COMM_WORLD, MPI_Info info=MPI_INFO_NULL)
Open an existing file.
NDArray is the datatype designed to hold any data, including scalars, multidimensional arrays,...