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
9NCTimeSeriesRiver::NCTimeSeriesRiver (const std::string a_file_name, const std::string a_field_name,
10 const std::string a_time_name,
11 const int a_nz, const int a_use_vert_integ) {
12 file_name = a_file_name;
13 time_name = a_time_name;
14 field_name = a_field_name;
15 nz = a_nz;
16 use_vert_integ = a_use_vert_integ;
17}
18
20 // open file
21 amrex::Print() << "Loading " << field_name << " from rivers NetCDF file " << file_name << std::endl;
22
23 // The time field can have any number of names, depending on the field.
24 // If not specified in input file (time_name.empty()) then set it by default
25 if (time_name.empty())
26 {
27 time_name = "river_time";
28 }
29
30 // Check units of time stamps; should be days
31 std::string unit_str = ReadNetCDFVarAttrStr(file_name, time_name, "units"); // works on proc 0
32 if (amrex::ParallelDescriptor::IOProcessor())
33 {
34 if (unit_str.find("days") == std::string::npos) {
35 amrex::Print() << "Units of river_time given as: " << unit_str << std::endl;
36 amrex::Abort("Units must be in days.");
37 }
38 }
39 // get times and put in array
40 using RARRAY = NDArray<amrex::Real>;
41 amrex::Vector<RARRAY> array_ts(1);
42 ReadNetCDFFile(file_name, {time_name}, array_ts); // filled only on proc 0
43 if (amrex::ParallelDescriptor::IOProcessor())
44 {
45 int ntimes_io = array_ts[0].get_vshape()[0];
46 for (int nt(0); nt < ntimes_io; nt++)
47 {
48 // Convert river time from days to seconds
49 river_times.push_back((*(array_ts[0].get_data() + nt)) * amrex::Real(60.0) * amrex::Real(60.0) * amrex::Real(24.0));
50 // amrex::Print() << "TIMES " << river_times[nt] << std::endl;
51 }
52 }
53 int ntimes = river_times.size();
54 amrex::Print() << "ntimes " << ntimes << std::endl;
55 int ioproc = amrex::ParallelDescriptor::IOProcessorNumber();
56 amrex::ParallelDescriptor::Bcast(&ntimes,1,ioproc);
57 if (!(amrex::ParallelDescriptor::IOProcessor())) {
58 river_times.resize(ntimes);
59 }
60 amrex::ParallelDescriptor::Bcast(river_times.data(), river_times.size(), ioproc);
61
62 auto ncf = ncutils::NCFile::open(file_name, NC_NOCLOBBER );
63 if (amrex::ParallelDescriptor::IOProcessor())
64 {
65 std::vector<MPI_Offset> shape = ncf.var(field_name).shape();
66 if (shape.size() == 2) {
67 has_z = 0;
68 nriv = shape[1];
69 } else if (shape.size() == 3) {
70 has_z = 1;
71 nriv = shape[2];
72 } else {
73 amrex::Abort("River field shape not 2 or 3");
74 }
75 }
76 amrex::ParallelDescriptor::Bcast(&has_z, 1, ioproc);
77 amrex::ParallelDescriptor::Bcast(&nriv, 1, ioproc);
78
79 amrex::Box vshape_box(amrex::IntVect(0,0,0),amrex::IntVect(nriv,0,nz));
80 if (!has_z && !use_vert_integ) {
81 amrex::Vector<amrex::FArrayBox*> NC_fabs;
82 amrex::Vector<std::string> NC_names;
83 amrex::Vector<enum NC_Data_Dims_Type> NC_dim_types;
84 amrex::Print() << "Reading in river_Vshape from " << file_name << std::endl;
85
86 fab_vshape = new amrex::FArrayBox();
87 NC_fabs.push_back(fab_vshape); NC_names.push_back("river_Vshape");
88 NC_dim_types.push_back(NC_Data_Dims_Type::BT_Riv);
89 BuildFABsFromNetCDFFile<amrex::FArrayBox,amrex::Real>(vshape_box, file_name, NC_names, NC_dim_types, NC_fabs);
90 }
91
92 nzbox = (use_vert_integ) ? 1 : nz;
93 amrex::Box riv_box(amrex::IntVect(0,0,0),amrex::IntVect(nriv,0,nzbox));
94#ifdef AMREX_USE_GPU
95 // It's possible there should be a different arena
96 fab_before = new amrex::FArrayBox(riv_box,1,amrex::The_Pinned_Arena());
97 fab_after = new amrex::FArrayBox(riv_box,1,amrex::The_Pinned_Arena());
98 fab_interp = new amrex::FArrayBox(riv_box,1,amrex::The_Pinned_Arena());
99#else
100 fab_before = new amrex::FArrayBox(riv_box,1);
101 fab_after = new amrex::FArrayBox(riv_box,1);
102 fab_interp = new amrex::FArrayBox(riv_box,1);
103#endif
104
105 // dummy initialization
106 i_time_before = -100;
107}
108
110 // Figure out time index:
111 AMREX_ASSERT(time >= river_times[0]);
112 AMREX_ASSERT(time <= river_times[river_times.size()-1]);
113 int i_time_before_old = i_time_before;
114 for (int nt=0; nt < river_times.size()-1; nt++) {
115 if ((river_times[nt] <= time) and (river_times[nt+1] >= time)) {
116 i_time_before = nt;
118 time_after = river_times[nt+1];
119 break;
120 }
121 }
122
123 int i_time_after = i_time_before + 1;
124 if (i_time_before_old + 1 == i_time_before) {
125 // swap data vectors so we only have to read in one MultiFab
126 std::swap(fab_before, fab_after);
127 read_in_at_time(fab_after, i_time_after);
128 } else if (i_time_before_old != i_time_before) {
129 read_in_at_time(fab_after, i_time_after);
131 }
132
133 amrex::Real dt = time_after - time_before;
134 amrex::Real time_before_copy = time_before;
135
136 amrex::Box fab_domain(amrex::IntVect(0,0,0), amrex::IntVect(nriv-1,0,nzbox-1));
137 auto interp_array = fab_interp->array();
138 auto before_array = fab_before->array();
139 auto after_array = fab_after->array();
140 amrex::ParallelFor(fab_domain, [=] AMREX_GPU_DEVICE (int r, int , int k) {
141 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;
142 });
143}
144
145void NCTimeSeriesRiver::read_in_at_time (amrex::FArrayBox* fab_dat, int itime) {
146 amrex::FArrayBox NC_fab;
147 amrex::Vector<amrex::FArrayBox*> NC_fabs;
148 amrex::Vector<std::string> NC_names;
149 amrex::Vector<enum NC_Data_Dims_Type> NC_dim_types;
150 // actual dims don't really matter here; only lower is used in call
151
152 amrex::Print() << "Reading in " << field_name << " at time index " << itime << " from " << file_name << std::endl;
153
154 NC_fabs.push_back(&NC_fab) ; NC_names.push_back(field_name);
155
156 if (has_z) {
157 NC_dim_types.push_back(NC_Data_Dims_Type::Time_BT_Riv);
158 } else {
159 NC_dim_types.push_back(NC_Data_Dims_Type::Time_Riv);
160 }
161
162 amrex::Box riv_domain(amrex::IntVect(0,0,0), amrex::IntVect(nriv-1,0,nz-1));
163 amrex::Box fab_domain(amrex::IntVect(0,0,0), amrex::IntVect(nriv-1,0,nzbox-1));
164
165 BuildFABsFromNetCDFFile<amrex::FArrayBox,amrex::Real>(riv_domain, file_name, NC_names, NC_dim_types, NC_fabs, true, itime);
166
167 auto dat_array = fab_dat->array();
168 auto tmp_array = NC_fabs[0]->const_array();
169 if (has_z || use_vert_integ) {
170 amrex::ParallelFor(fab_domain, [=] AMREX_GPU_DEVICE (int r, int , int k) {
171 dat_array(r,0,k) = tmp_array(r,0,k);
172 });
173 } else {
174 auto array_vshape = fab_vshape->const_array();
175 amrex::ParallelFor(fab_domain, [=] AMREX_GPU_DEVICE (int r, int , int k) {
176 dat_array(r,0,k) = tmp_array(r,0,0) * array_vshape(r,0,k);
177 });
178 }
179}
180#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.
std::string file_name
File name to read from.
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.
std::string time_name
Field name for time series in netcdf file.
amrex::Vector< amrex::Real > river_times
Time points in netcdf file.
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::FArrayBox * fab_interp
Container for interpolated data; Only used if save_interpolated == true.
NCTimeSeriesRiver(const std::string a_file_name, const std::string a_field_name, const std::string a_time_name, const int a_nz, const int a_use_vert_integ=0)
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,...