REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_NCTimeSeries.cpp
Go to the documentation of this file.
2#include "REMORA_NCFile.H"
3
4#include "AMReX_FillPatchUtil.H"
5#include "AMReX_Interpolater.H"
6#include "AMReX_ParallelDescriptor.H"
7
8#include <string>
9
10#ifdef REMORA_USE_NETCDF
11/**
12 * @param[in ] a_file_name file name to read from
13 * @param[in ] a_field_name name of field to read in
14 * @param[in ] a_time_name name of time variable in NetCDF file
15 * @param[in ] a_domain simulation domain
16 * @param[inout] a_mf_var MultiFab of data to either store into or reference for dimensions
17 * @param[in ] a_is2d Whether the variable we're working with is 2D
18 * @param[in ] a_save_interpolated Whether the interpolated value should be saved internally
19 */
20NCTimeSeries::NCTimeSeries (const std::string a_file_name, const std::string a_field_name,
21 const std::string a_time_name,
22 const amrex::Box& a_domain,
23 amrex::MultiFab* a_mf_var, bool a_is2d, bool a_save_interpolated) {
24 file_name = a_file_name;
25 time_name = a_time_name;
26 field_name = a_field_name;
27 domain = a_domain;
28 mf_var = a_mf_var;
29 is2d = a_is2d;
30 save_interpolated = a_save_interpolated;
31}
32
34 // open file
35 amrex::Print() << "Loading " << field_name << " from NetCDF file " << file_name << std::endl;
36
37 // The time field can have any number of names, depending on the field.
38 // If not specified in input file (time_name.empty()) then set it by default
39 if (time_name.empty())
40 {
41 if (field_name.find("wind") != std::string::npos) {
42 time_name = "wind_time";
43 } else if ((field_name.find("str") != std::string::npos) and (field_name[0] == 's')) {
44 time_name = "sms_time";
45 } else if ((field_name.find("str") != std::string::npos) and (field_name[0] == 'b')) {
46 time_name = "bms_time";
47 } else {
48 time_name = "ocean_time";
49 }
50 }
51
52 // Check units of time stamps; should be days
53 std::string unit_str = ReadNetCDFVarAttrStr(file_name, time_name, "units"); // works on proc 0
54 if (amrex::ParallelDescriptor::IOProcessor())
55 {
56 if (unit_str.find("days") == std::string::npos) {
57 amrex::Print() << "Units of ocean_time given as: " << unit_str << std::endl;
58 amrex::Abort("Units must be in days.");
59 }
60 }
61 // get times and put in array
62 using RARRAY = NDArray<amrex::Real>;
63 amrex::Vector<RARRAY> array_ts(1);
64 ReadNetCDFFile(file_name, {time_name}, array_ts); // filled only on proc 0
65 if (amrex::ParallelDescriptor::IOProcessor())
66 {
67 int ntimes_io = array_ts[0].get_vshape()[0];
68 for (int nt(0); nt < ntimes_io; nt++)
69 {
70 // Convert ocean time from days to seconds
71 ocean_times.push_back((*(array_ts[0].get_data() + nt)) * amrex::Real(60.0) * amrex::Real(60.0) * amrex::Real(24.0));
72 // amrex::Print() << "TIMES " << ocean_times[nt] << std::endl;
73 }
74 }
75 int ntimes = ocean_times.size();
76 // Only do checks on IO processor since ocean_times isn't populated on other ranks yet
77 if (amrex::ParallelDescriptor::IOProcessor()) {
78 AMREX_ASSERT(std::is_sorted(ocean_times.begin(), ocean_times.end()));
79 if (ntimes <= 1) {
80 amrex::Error("Time series data must be given at at least two times");
81 }
82 }
83 int ioproc = amrex::ParallelDescriptor::IOProcessorNumber();
84 amrex::ParallelDescriptor::Bcast(&ntimes,1,ioproc);
85 if (!(amrex::ParallelDescriptor::IOProcessor())) {
86 ocean_times.resize(ntimes);
87 }
88 amrex::ParallelDescriptor::Bcast(ocean_times.data(), ocean_times.size(), ioproc);
89
90 // Initialize MultiFabs
91 // NetCDF data is always read and temporally interpolated on level 0.
92 mf_before = new amrex::MultiFab(mf_var->boxArray(), mf_var->DistributionMap(), 1, mf_var->nGrowVect());
93 mf_after = new amrex::MultiFab(mf_var->boxArray(), mf_var->DistributionMap(), 1, mf_var->nGrowVect());
94 mf_interp_lev0 = new amrex::MultiFab(mf_var->boxArray(), mf_var->DistributionMap(), 1, mf_var->nGrowVect());
95
96 // dummy initialization
97 i_time_before = -100;
98}
99
100/**
101 * @param time time to interpolate to
102 */
103void NCTimeSeries::update_interpolated_to_time (amrex::Real time, int lev,
104 amrex::MultiFab* mf_lev,
105 const amrex::Vector<amrex::Geometry>& geom,
106 const amrex::Vector<amrex::IntVect>& ref_ratio) {
107 // Figure out time index:
108 AMREX_ASSERT(time >= ocean_times[0]);
109 AMREX_ASSERT(time <= ocean_times[ocean_times.size()-1]);
110 int i_time_before_old = i_time_before;
111 for (int nt=0; nt < ocean_times.size()-1; nt++) {
112 if ((ocean_times[nt] <= time) and (ocean_times[nt+1] >= time)) {
113 i_time_before = nt;
115 time_after = ocean_times[nt+1];
116 break;
117 }
118 }
119
120 int i_time_after = i_time_before + 1;
121 if (i_time_before_old + 1 == i_time_before) {
122 // swap multifabs so we only have to read in one MultiFab
123 std::swap(mf_before, mf_after);
124 read_in_at_time(mf_after, i_time_after);
125 } else if (i_time_before_old != i_time_before) {
126 read_in_at_time(mf_after, i_time_after);
128 }
129
130 amrex::Real dt = time_after - time_before;
131
132 auto nodality = mf_interp_lev0->ixType();
133
134#ifdef AMREX_USE_OMP
135#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
136#endif
137 for (amrex::MFIter mfi(*mf_interp_lev0,true); mfi.isValid(); ++mfi) {
138 // Adjust box to match ROMS grid
139 amrex::Box bx = mfi.growntilebox(amrex::IntVect(1-nodality[0],1-nodality[1],0));
140
141 amrex::Real time_before_copy = time_before;
142
143 // Temporal interpolation is done once on level 0.
144 amrex::MultiFab* mf_to_fill = mf_interp_lev0;
145 amrex::Array4<amrex::Real> to_fill = mf_to_fill->array(mfi);
146 amrex::Array4<const amrex::Real> before = mf_before->const_array(mfi);
147 amrex::Array4<const amrex::Real> after = mf_after->const_array(mfi);
148 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
149 {
150 to_fill(i,j,k) = before(i,j,k) + (time - time_before_copy) * (after(i,j,k) - before(i,j,k)) / dt;
151 });
152 }
153
154 amrex::MultiFab* mf_to_fill_lev = mf_lev;
155 if (save_interpolated) {
156 if (mf_interpolated_lev.size() <= static_cast<std::size_t>(lev)) {
157 mf_interpolated_lev.resize(lev+1);
158 }
159 if (!mf_interpolated_lev[lev] ||
160 mf_interpolated_lev[lev]->boxArray() != mf_lev->boxArray() ||
161 mf_interpolated_lev[lev]->nGrowVect() != mf_lev->nGrowVect()) {
162 mf_interpolated_lev[lev] = std::make_unique<amrex::MultiFab>(
163 mf_lev->boxArray(), mf_lev->DistributionMap(), 1, mf_lev->nGrowVect());
164 }
165 mf_to_fill_lev = mf_interpolated_lev[lev].get();
166 }
167
168 if (lev == 0) {
169 amrex::MultiFab::Copy(*mf_to_fill_lev, *mf_interp_lev0, 0, 0, 1, mf_to_fill_lev->nGrowVect());
170 return;
171 }
172
173 amrex::PhysBCFunctNoOp null_bc_for_fill;
174 amrex::Vector<amrex::BCRec> null_dom_bcs(1);
175 for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
176 null_dom_bcs[0].setLo(dir, amrex::BCType::int_dir);
177 null_dom_bcs[0].setHi(dir, amrex::BCType::int_dir);
178 }
179
180 amrex::Interpolater* mapper = nullptr;
181 const auto& idx_type = mf_to_fill_lev->ixType();
182 if (idx_type == amrex::IndexType(amrex::IntVect(0,0,0))) {
183 mapper = &amrex::cell_cons_interp;
184 } else {
185 mapper = &amrex::face_cons_linear_interp;
186 }
187
188 amrex::InterpFromCoarseLevel(*mf_to_fill_lev, amrex::Real(0.0), *mf_interp_lev0,
189 0, 0, 1,
190 geom[0], geom[lev],
192 cumulative_ref_ratio(lev, ref_ratio),
193 mapper, null_dom_bcs, 0);
194
195 mf_to_fill_lev->FillBoundary(geom[lev].periodicity());
196}
197
198amrex::IntVect
200 const amrex::Vector<amrex::IntVect>& ref_ratio) const {
201 amrex::IntVect rr(1,1,1);
202 for (int l = 0; l < lev; ++l) {
203 rr[0] *= ref_ratio[l][0];
204 rr[1] *= ref_ratio[l][1];
205 rr[2] *= ref_ratio[l][2];
206 }
207 return rr;
208}
209
210const amrex::MultiFab*
212 AMREX_ALWAYS_ASSERT(save_interpolated);
213 AMREX_ALWAYS_ASSERT(lev < static_cast<int>(mf_interpolated_lev.size()));
214 AMREX_ALWAYS_ASSERT(mf_interpolated_lev[lev]);
215 return mf_interpolated_lev[lev].get();
216}
217
218/**
219 * @param[inout] mf multifab to store time step data into
220 * @param[in ] itime index of time step to read from file
221 */
222void NCTimeSeries::read_in_at_time (amrex::MultiFab* mf, int itime) {
223 // This all assumes that we're on level 0 with only one boxes_at_level
224 amrex::FArrayBox NC_fab;
225 amrex::Vector<amrex::FArrayBox*> NC_fabs;
226 amrex::Vector<std::string> NC_names;
227 amrex::Vector<enum NC_Data_Dims_Type> NC_dim_types;
228
229 amrex::Print() << "Reading in " << field_name << " at time index " << itime << " from " << file_name << std::endl;
230
231 NC_fabs.push_back(&NC_fab) ; NC_names.push_back(field_name);
232
233 if (is2d) {
234 NC_dim_types.push_back(NC_Data_Dims_Type::Time_SN_WE);
235 } else {
236 NC_dim_types.push_back(NC_Data_Dims_Type::Time_BT_SN_WE);
237 }
238
239 BuildFABsFromNetCDFFile<amrex::FArrayBox,amrex::Real>(domain, file_name, NC_names, NC_dim_types, NC_fabs, true, itime);
240
241#ifdef _OPENMP
242#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
243#endif
244 {
245 // Don't tile this since we are operating on full FABs in this routine
246 for ( amrex::MFIter mfi(*mf, false); mfi.isValid(); ++mfi )
247 {
248 amrex::FArrayBox &fab = (*mf)[mfi];
249
250 //
251 // FArrayBox to FArrayBox copy does "copy on intersection"
252 // This only works here because we have broadcast the FArrayBox of data from the netcdf file to all ranks
253
254 fab.template copy<amrex::RunOn::Device>(NC_fab);
255 } // mf
256 } // omp
257}
258#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.
PhysBCFunctNoOp null_bc_for_fill
amrex::Vector< std::unique_ptr< amrex::MultiFab > > mf_interpolated_lev
Interpolated data on each requested AMR level if save_interpolated=true.
bool is2d
Whether the field we're reading in is 2d.
void update_interpolated_to_time(amrex::Real time, int lev, amrex::MultiFab *mf_lev, const amrex::Vector< amrex::Geometry > &geom, const amrex::Vector< amrex::IntVect > &ref_ratio)
Calculate interpolated values at time and fill data for level lev.
void read_in_at_time(amrex::MultiFab *mf, int itime)
Read in data from file at time index itime and fill into mf.
amrex::MultiFab * mf_interp_lev0
Multifab storing temporally interpolated data on level 0.
amrex::Real time_after
Time in ocean_times immediately after the last time interpolated to.
amrex::Real time_before
Time in ocean_times immediately before the last time interpolated to.
amrex::Box domain
Domain.
void Initialize()
Read in time array from file and allocate data arrays.
amrex::MultiFab * mf_before
Multifab to store data at time_before.
int i_time_before
Time index immediately before the last time interpolated to.
std::string field_name
Field name in netcdf file.
amrex::IntVect cumulative_ref_ratio(int lev, const amrex::Vector< amrex::IntVect > &ref_ratio) const
Build cumulative refinement ratio from level 0 to lev.
amrex::MultiFab * mf_after
Multifab to store data at time_after.
std::string time_name
Field name for time series in netcdf file.
std::string file_name
File name to read from.
const amrex::MultiFab * get_interpolated_mf(int lev) const
Access interpolated data saved for a specific level.
amrex::Vector< amrex::Real > ocean_times
Time points in netcdf file.
amrex::MultiFab * mf_var
NCTimeSeries(const std::string a_file_name, const std::string a_field_name, const std::string a_time_name, const amrex::Box &a_domain, amrex::MultiFab *a_mf_var, bool a_is2d, bool a_save_interpolated)
Constructor.
NDArray is the datatype designed to hold any data, including scalars, multidimensional arrays,...