REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_NCTimeSeriesBoundary.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 * @param[in ] a_lev level at which we will store the data
11 * @param[in ] a_geom Vector of Geometry objects at all levels
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_index_type nodality of data field
16 * @param[in ] a_var_need_data array over boundaries of flags that indicate whether we need data for the variable
17 * @param[in ] a_is2d Whether the variable we're working with is 2D
18 * @param[in ] a_rx refinement ratio in x relative to level 0
19 * @param[in ] a_ry refinement ratio in y relative to level 0
20 */
21NCTimeSeriesBoundary::NCTimeSeriesBoundary (int a_lev, const amrex::Vector<amrex::Geometry> a_geom,
22 const amrex::Vector<std::string>& a_file_names, const std::string a_field_name,
23 const std::string a_time_name,
24 const amrex::IntVect a_index_type,
25 const amrex::GpuArray<bool, AMREX_SPACEDIM*2>* a_var_need_data,
26 bool a_is2d, int a_rx, int a_ry)
27{
28 m_lev = a_lev;
29 m_geom = a_geom;
30 file_names.assign(a_file_names.begin(), a_file_names.end());
31 time_name = a_time_name;
32 field_name = a_field_name;
33 domain = a_geom[a_lev].Domain();
34 index_type = a_index_type;
35 var_need_data = *a_var_need_data;
36 is2d = a_is2d;
37 m_rx = a_rx;
38 m_ry = a_ry;
39}
40
42{
43 // Initialize Fabs
44 amrex::Arena* Arena_Used = amrex::The_Arena();
45#ifdef AMREX_USE_GPU
46 Arena_Used = amrex::The_Pinned_Arena();
47#endif
48
49 // open file
50 amrex::Print() << "Setting up boundary data for " << field_name << " coming from NetCDF file " << std::endl;
51
52 // The time field can have any number of names, depending on the field.
53 // If not specified in input file (time_name.empty()) then set it by default
54 if (time_name.empty())
55 {
56 time_name = "ocean_time";
57 }
58
59 for (int ifile = 0; ifile < file_names.size(); ifile++) {
60 std::string file_name = file_names[ifile];
61 // Check units of time stamps; should be days
62 std::string unit_str = ReadNetCDFVarAttrStr(file_name, time_name, "units"); // works on proc 0
63 if (amrex::ParallelDescriptor::IOProcessor())
64 {
65 if (unit_str.find("days") == std::string::npos) {
66 amrex::Print() << "Units of ocean_time given as: " << unit_str << std::endl;
67 amrex::Abort("Units must be in days.");
68 }
69 }
70 // get times and put in array
71 using RARRAY = NDArray<amrex::Real>;
72 amrex::Vector<RARRAY> array_ts(1);
73 ReadNetCDFFile(file_name, {time_name}, array_ts); // filled only on proc 0
74 if (amrex::ParallelDescriptor::IOProcessor())
75 {
76 int ntimes_io = array_ts[0].get_vshape()[0];
77 for (int nt(0); nt < ntimes_io; nt++)
78 {
79 // Convert ocean time from days to seconds
80 bry_times.push_back((*(array_ts[0].get_data() + nt)) * amrex::Real(60.0) * amrex::Real(60.0) * amrex::Real(24.0));
81 file_for_time.push_back(ifile);
82 file_itime_offset.push_back(nt);
83 // amrex::Print() << "TIMES " << bry_times[nt] << std::endl;
84 }
85 }
86 }
87
88 int ntimes = bry_times.size();
89 // Only do checks on IO processor since bry_times isn't populated on other ranks yet
90 if (amrex::ParallelDescriptor::IOProcessor()) {
91 AMREX_ASSERT(std::is_sorted(bry_times.begin(), bry_times.end()));
92 if (ntimes <= 1) {
93 amrex::Error("Time series of boundary data must be given at at least two times");
94 }
95 }
96 int ioproc = amrex::ParallelDescriptor::IOProcessorNumber();
97 amrex::ParallelDescriptor::Bcast(&ntimes,1,ioproc);
98 if (!(amrex::ParallelDescriptor::IOProcessor())) {
99 bry_times.resize(ntimes);
100 file_for_time.resize(ntimes);
101 file_itime_offset.resize(ntimes);
102 }
103 amrex::ParallelDescriptor::Bcast(bry_times.data(), bry_times.size(), ioproc);
104 amrex::ParallelDescriptor::Bcast(file_for_time.data(), file_for_time.size(), ioproc);
105 amrex::ParallelDescriptor::Bcast(file_itime_offset.data(), file_itime_offset.size(), ioproc);
106
107 const auto& lo = domain.loVect();
108 const auto& hi = domain.hiVect();
109
110 amrex::Box xlo_bx(amrex::IntVect(lo[0]+index_type[0]-1, lo[1]+index_type[1]-1, lo[2]),
111 amrex::IntVect(lo[0]+index_type[0]-1, hi[1]+1 , hi[2]), index_type);
112 amrex::Box xhi_bx(amrex::IntVect(hi[0]+1 , lo[1]+index_type[1]-1, lo[2]),
113 amrex::IntVect(hi[0]+1 , hi[1]+1 , hi[2]), index_type);
114 amrex::Box ylo_bx(amrex::IntVect(lo[0]+index_type[0]-1, lo[1]+index_type[1]-1, lo[2]),
115 amrex::IntVect(hi[0]+1 , lo[1]+index_type[1]-1, hi[2]), index_type);
116 amrex::Box yhi_bx(amrex::IntVect(lo[0]+index_type[0]-1, hi[1]+1 , lo[2]),
117 amrex::IntVect(hi[0]+1 , hi[1]+1 , hi[2]), index_type);
118 if (is2d) {
119 xlo_bx.makeSlab(2,0);
120 xhi_bx.makeSlab(2,0);
121 ylo_bx.makeSlab(2,0);
122 yhi_bx.makeSlab(2,0);
123 }
124
125 // amrex::Print() << xlo_bx << " " << xhi_bx << " " << ylo_bx << " " << yhi_bx << std::endl;
126
127 xlo_dat_before = amrex::FArrayBox(xlo_bx, 1, Arena_Used);
128 xhi_dat_before = amrex::FArrayBox(xhi_bx, 1, Arena_Used);
129 ylo_dat_before = amrex::FArrayBox(ylo_bx, 1, Arena_Used);
130 yhi_dat_before = amrex::FArrayBox(yhi_bx, 1, Arena_Used);
131
132 xlo_dat_after = amrex::FArrayBox(xlo_bx, 1, Arena_Used);
133 xhi_dat_after = amrex::FArrayBox(xhi_bx, 1, Arena_Used);
134 ylo_dat_after = amrex::FArrayBox(ylo_bx, 1, Arena_Used);
135 yhi_dat_after = amrex::FArrayBox(yhi_bx, 1, Arena_Used);
136
137 xlo_dat_interp = amrex::FArrayBox(xlo_bx, 1, Arena_Used);
138 xhi_dat_interp = amrex::FArrayBox(xhi_bx, 1, Arena_Used);
139 ylo_dat_interp = amrex::FArrayBox(ylo_bx, 1, Arena_Used);
140 yhi_dat_interp = amrex::FArrayBox(yhi_bx, 1, Arena_Used);
141
142 // dummy initialization
143 i_time_before = -100;
144
145 if (var_need_data[amrex::Orientation(amrex::Direction::x,amrex::Orientation::low)] == true) {
146 nc_var_names.push_back(field_name + "_west");
147 }
148 if (var_need_data[amrex::Orientation(amrex::Direction::x,amrex::Orientation::high)] == true) {
149 nc_var_names.push_back(field_name + "_east");
150 }
151 if (var_need_data[amrex::Orientation(amrex::Direction::y,amrex::Orientation::low)] == true) {
152 nc_var_names.push_back(field_name + "_south");
153 }
154 if (var_need_data[amrex::Orientation(amrex::Direction::y,amrex::Orientation::high)] == true) {
155 nc_var_names.push_back(field_name + "_north");
156 }
157}
158
159/**
160 * @param time time to interpolate to
161 */
163{
164 // Initialize Fabs
165 amrex::Arena* Arena_Used = amrex::The_Arena();
166#ifdef AMREX_USE_GPU
167 Arena_Used = amrex::The_Pinned_Arena();
168#endif
169
170 // Figure out time index:
171 AMREX_ASSERT(time >= bry_times[0]);
172 AMREX_ASSERT(time <= bry_times[bry_times.size()-1]);
173 int i_time_before_old = i_time_before;
174 for (int nt=0; nt < bry_times.size()-1; nt++) {
175 if ((bry_times[nt] <= time) and (bry_times[nt+1] >= time)) {
176 i_time_before = nt;
178 time_after = bry_times[nt+1];
179 break;
180 }
181 }
182
183 int i_time_after = i_time_before + 1;
184
185 amrex::FArrayBox crse_xlo_dat;
186 amrex::FArrayBox crse_xhi_dat;
187 amrex::FArrayBox crse_ylo_dat;
188 amrex::FArrayBox crse_yhi_dat;
189
190 if (m_lev > 0) {
191 const auto& crse_lo = m_geom[0].Domain().loVect();
192 const auto& crse_hi = m_geom[0].Domain().hiVect();
193
194 amrex::Box crse_xlo_bx(amrex::IntVect(crse_lo[0]+index_type[0]-1, crse_lo[1]+index_type[1]-1, crse_lo[2]),
195 amrex::IntVect(crse_lo[0]+index_type[0]-1, crse_hi[1]+1 , crse_hi[2]), index_type);
196 amrex::Box crse_xhi_bx(amrex::IntVect(crse_hi[0]+1 , crse_lo[1]+index_type[1]-1, crse_lo[2]),
197 amrex::IntVect(crse_hi[0]+1 , crse_hi[1]+1 , crse_hi[2]), index_type);
198 amrex::Box crse_ylo_bx(amrex::IntVect(crse_lo[0]+index_type[0]-1, crse_lo[1]+index_type[1]-1, crse_lo[2]),
199 amrex::IntVect(crse_hi[0]+1 , crse_lo[1]+index_type[1]-1, crse_hi[2]), index_type);
200 amrex::Box crse_yhi_bx(amrex::IntVect(crse_lo[0]+index_type[0]-1, crse_hi[1]+1 , crse_lo[2]),
201 amrex::IntVect(crse_hi[0]+1 , crse_hi[1]+1 , crse_hi[2]), index_type);
202 if (is2d) {
203 crse_xlo_bx.makeSlab(2,0);
204 crse_xhi_bx.makeSlab(2,0);
205 crse_ylo_bx.makeSlab(2,0);
206 crse_yhi_bx.makeSlab(2,0);
207 }
208
209 crse_xlo_dat = amrex::FArrayBox(crse_xlo_bx, 1, Arena_Used);
210 crse_xhi_dat = amrex::FArrayBox(crse_xhi_bx, 1, Arena_Used);
211 crse_ylo_dat = amrex::FArrayBox(crse_ylo_bx, 1, Arena_Used);
212 crse_yhi_dat = amrex::FArrayBox(crse_yhi_bx, 1, Arena_Used);
213 }
214
215 if (i_time_before_old + 1 == i_time_before) {
216 // swap multifabs so we only have to read in one MultiFab
217 std::swap(xlo_dat_before, xlo_dat_after);
218 std::swap(xhi_dat_before, xhi_dat_after);
219 std::swap(ylo_dat_before, ylo_dat_after);
220 std::swap(yhi_dat_before, yhi_dat_after);
221
222 amrex::Print() << "Reading in " << field_name << " at (after) time " << i_time_after << std::endl;
223
224 if (m_lev == 0) {
226 } else {
227
228 read_in_at_time(crse_xlo_dat, crse_xhi_dat, crse_ylo_dat, crse_yhi_dat, i_time_after);
229
230 interp_fab(crse_xlo_dat, xlo_dat_after);
231 interp_fab(crse_xhi_dat, xhi_dat_after);
232 interp_fab(crse_ylo_dat, ylo_dat_after);
233 interp_fab(crse_yhi_dat, yhi_dat_after);
234
235 }
236
237 } else if (i_time_before_old != i_time_before) {
238
239
240 if (m_lev == 0) {
241 amrex::Print() << "Reading in " << field_name << " at (before) time " << i_time_before << std::endl;
243 amrex::Print() << "Reading in " << field_name << " at (after) time " << i_time_after << std::endl;
245 } else {
246
247 amrex::Print() << "Reading in " << field_name << " at (before) time " << i_time_before << std::endl;
248
249 read_in_at_time(crse_xlo_dat, crse_xhi_dat, crse_ylo_dat,crse_yhi_dat, i_time_before);
250
251 interp_fab(crse_xlo_dat, xlo_dat_before);
252 interp_fab(crse_xhi_dat, xhi_dat_before);
253 interp_fab(crse_ylo_dat, ylo_dat_before);
254 interp_fab(crse_yhi_dat, yhi_dat_before);
255
256 amrex::Print() << "Reading in " << field_name << " at time " << i_time_after << std::endl;
257
258 read_in_at_time(crse_xlo_dat, crse_xhi_dat, crse_ylo_dat, crse_yhi_dat, i_time_after);
259
260 interp_fab(crse_xlo_dat, xlo_dat_after);
261 interp_fab(crse_xhi_dat, xhi_dat_after);
262 interp_fab(crse_ylo_dat, ylo_dat_after);
263 interp_fab(crse_yhi_dat, yhi_dat_after);
264 } // lev
265 } // i_time
266
267 amrex::Real dt = time_after - time_before;
268 amrex::Real time_before_copy = time_before;
269
270 amrex::Array4<amrex::Real> xlo_interp_arr = xlo_dat_interp.array();
271 amrex::Array4<amrex::Real> xhi_interp_arr = xhi_dat_interp.array();
272 amrex::Array4<amrex::Real> ylo_interp_arr = ylo_dat_interp.array();
273 amrex::Array4<amrex::Real> yhi_interp_arr = yhi_dat_interp.array();
274
275 amrex::Array4<amrex::Real> xlo_before_arr = xlo_dat_before.array();
276 amrex::Array4<amrex::Real> xhi_before_arr = xhi_dat_before.array();
277 amrex::Array4<amrex::Real> ylo_before_arr = ylo_dat_before.array();
278 amrex::Array4<amrex::Real> yhi_before_arr = yhi_dat_before.array();
279
280 amrex::Array4<amrex::Real> xlo_after_arr = xlo_dat_after.array();
281 amrex::Array4<amrex::Real> xhi_after_arr = xhi_dat_after.array();
282 amrex::Array4<amrex::Real> ylo_after_arr = ylo_dat_after.array();
283 amrex::Array4<amrex::Real> yhi_after_arr = yhi_dat_after.array();
284
285 if (var_need_data[amrex::Orientation(amrex::Direction::x,amrex::Orientation::low)] == true) {
286 amrex::ParallelFor(xlo_dat_interp.box(), [=] AMREX_GPU_DEVICE (int i, int j, int k)
287 {
288 xlo_interp_arr(i,j,k) = xlo_before_arr(i,j,k) + (time - time_before_copy) * (xlo_after_arr(i,j,k) - xlo_before_arr(i,j,k)) / dt;
289 });
290 }
291
292 if (var_need_data[amrex::Orientation(amrex::Direction::x,amrex::Orientation::high)] == true) {
293 amrex::ParallelFor(xhi_dat_interp.box(), [=] AMREX_GPU_DEVICE (int i, int j, int k)
294 {
295 xhi_interp_arr(i,j,k) = xhi_before_arr(i,j,k) + (time - time_before_copy) * (xhi_after_arr(i,j,k) - xhi_before_arr(i,j,k)) / dt;
296 });
297 }
298
299 if (var_need_data[amrex::Orientation(amrex::Direction::y,amrex::Orientation::low)] == true) {
300 amrex::ParallelFor(ylo_dat_interp.box(), [=] AMREX_GPU_DEVICE (int i, int j, int k)
301 {
302 ylo_interp_arr(i,j,k) = ylo_before_arr(i,j,k) + (time - time_before_copy) * (ylo_after_arr(i,j,k) - ylo_before_arr(i,j,k)) / dt;
303 });
304 }
305
306 if (var_need_data[amrex::Orientation(amrex::Direction::y,amrex::Orientation::high)] == true) {
307 amrex::ParallelFor(yhi_dat_interp.box(), [=] AMREX_GPU_DEVICE (int i, int j, int k)
308 {
309 yhi_interp_arr(i,j,k) = yhi_before_arr(i,j,k) + (time - time_before_copy) * (yhi_after_arr(i,j,k) - yhi_before_arr(i,j,k)) / dt;
310 });
311 }
312}
313
314/**
315 * @param[inout] fab_xlo fab to store xlo boundary data into
316 * @param[inout] fab_xhi fab to store xhi boundary data into
317 * @param[inout] fab_ylo fab to store ylo boundary data into
318 * @param[inout] fab_yhi fab to store yhi boundary data into
319 * @param[in ] itime index of time step to read from file
320 */
321void NCTimeSeriesBoundary::read_in_at_time (amrex::FArrayBox& fab_xlo,
322 amrex::FArrayBox& fab_xhi,
323 amrex::FArrayBox& fab_ylo,
324 amrex::FArrayBox& fab_yhi,
325 int itime) {
326 using RARRAY = NDArray<amrex::Real>;
327 amrex::Vector<RARRAY> arrays(nc_var_names.size());
328
329 amrex::Print() << "Actually reading in " << field_name << " at time " << bry_times[itime] << std::endl;
330 std::string nc_bdry_file = file_names[file_for_time[itime]];
331 int itime_offset = file_itime_offset[itime];
332 ReadNetCDFFile(nc_bdry_file, nc_var_names, arrays, true, itime_offset); // does work on proc 0 only
333
334 for (int iv=0; iv < nc_var_names.size(); iv++) {
335 if (amrex::ParallelDescriptor::IOProcessor())
336 {
337 std::string last4 = nc_var_names[iv].substr(nc_var_names[iv].size()-4, 4);
338 std::string last5 = nc_var_names[iv].substr(nc_var_names[iv].size()-5, 5);
339 int nx, ny, nz, n_plane;
340 int i, j, k, ioff, joff;
341 if (last4 == "west") {
342 amrex::Box my_box = fab_xlo.box();
343 if (is2d) {
344 nz = 1;
345 ny = arrays[iv].get_vshape()[1];
346 } else {
347 nz = arrays[iv].get_vshape()[1];
348 ny = arrays[iv].get_vshape()[2];
349 }
350 n_plane = ny * nz;
351 AMREX_ALWAYS_ASSERT(my_box.numPts() == n_plane);
352
353 i = my_box.smallEnd()[0];
354 joff = my_box.smallEnd()[1];
355
356 amrex::Array4<amrex::Real> fab_arr = fab_xlo.array();
357
358 for (int n(0); n < n_plane; n++) {
359 k = n / ny;
360 j = n - (k * ny);
361 fab_arr(i, j+joff, k, 0) = static_cast<amrex::Real>(*(arrays[iv].get_data() + n));
362 }
363 } else if (last4 == "east") {
364 amrex::Box my_box = fab_xhi.box();
365 if (is2d) {
366 nz = 1;
367 ny = arrays[iv].get_vshape()[1];
368 } else {
369 nz = arrays[iv].get_vshape()[1];
370 ny = arrays[iv].get_vshape()[2];
371 }
372 n_plane = ny * nz;
373 AMREX_ALWAYS_ASSERT(my_box.numPts() == n_plane);
374
375 i = my_box.smallEnd()[0];
376 joff = my_box.smallEnd()[1];
377
378 amrex::Array4<amrex::Real> fab_arr = fab_xhi.array();
379 for (int n(0); n < n_plane; n++) {
380 k = n / ny;
381 j = n - (k * ny);
382 fab_arr(i, j+joff, k, 0) = static_cast<amrex::Real>(*(arrays[iv].get_data() + n));
383 }
384 } else if (last5 == "south") {
385 amrex::Box my_box = fab_ylo.box();
386 if (is2d) {
387 nz = 1;
388 nx = arrays[iv].get_vshape()[1];
389 } else {
390 nz = arrays[iv].get_vshape()[1];
391 nx = arrays[iv].get_vshape()[2];
392 }
393 n_plane = nx * nz;
394 AMREX_ALWAYS_ASSERT(my_box.numPts() == n_plane);
395
396 j = my_box.smallEnd()[1];
397 ioff = my_box.smallEnd()[0];
398
399 amrex::Array4<amrex::Real> fab_arr = fab_ylo.array();
400 for (int n(0); n < n_plane; n++) {
401 k = n / nx;
402 i = n - (k * nx);
403 fab_arr(i+ioff, j, k, 0) = static_cast<amrex::Real>(*(arrays[iv].get_data() + n));
404 }
405 } else if (last5 == "north") {
406 amrex::Box my_box = fab_yhi.box();
407 if (is2d) {
408 nz = 1;
409 nx = arrays[iv].get_vshape()[1];
410 } else {
411 nz = arrays[iv].get_vshape()[1];
412 nx = arrays[iv].get_vshape()[2];
413 }
414 n_plane = nx * nz;
415 AMREX_ALWAYS_ASSERT(my_box.numPts() == n_plane);
416
417 j = my_box.smallEnd()[1];
418 ioff = my_box.smallEnd()[0];
419
420 amrex::Array4<amrex::Real> fab_arr = fab_yhi.array();
421 for (int n(0); n < n_plane; n++) {
422 k = n / nx;
423 i = n - (k * nx);
424 fab_arr(i+ioff, j, k, 0) = static_cast<amrex::Real>(*(arrays[iv].get_data() + n));
425 }
426 }
427 }
428 }
429
430 amrex::ParallelDescriptor::Barrier();
431
432 // When an FArrayBox is built, space is allocated on every rank. However, we only
433 // filled the data in these FABs on the IOProcessor. So here we broadcast
434 // the data to every rank.
435 int ioproc = amrex::ParallelDescriptor::IOProcessorNumber(); // I/O rank
436 amrex::ParallelDescriptor::Bcast(fab_xlo.dataPtr(),fab_xlo.box().numPts(),ioproc);
437 amrex::ParallelDescriptor::Bcast(fab_xhi.dataPtr(),fab_xhi.box().numPts(),ioproc);
438 amrex::ParallelDescriptor::Bcast(fab_ylo.dataPtr(),fab_ylo.box().numPts(),ioproc);
439 amrex::ParallelDescriptor::Bcast(fab_yhi.dataPtr(),fab_yhi.box().numPts(),ioproc);
440
441 amrex::Print() << "DONE reading in " << field_name << " at time " << bry_times[itime] << std::endl;
442}
443
444/**
445 * @param[in ] dat_crse fab of coarse data to interpolate from
446 * @param[ out] dat_fine fab of fine data to interpoalte to
447 */
448void NCTimeSeriesBoundary::interp_fab(amrex::FArrayBox& dat_crse, amrex::FArrayBox& dat_fine)
449{
450 amrex::Array4<amrex::Real> crse_arr = dat_crse.array();
451 amrex::Array4<amrex::Real> fine_arr = dat_fine.array();
452
453 const auto& bhi = ubound(dat_crse.box());
454
455 int rx = m_rx;
456 int ry = m_ry;
457
458 amrex::Real xfac = 1.0 / static_cast<amrex::Real>(rx);
459 amrex::Real yfac = 1.0 / static_cast<amrex::Real>(ry);
460
461 // Doing box on x-face
462 if (dat_crse.box().length(0) == 1) {
463 // amrex::Print() << "DOING INTERP ON XFACE " << dat_crse.box() << " " << dat_fine.box() << std::endl;
464 if (dat_crse.box().ixType()[1] == 0) {
465 amrex::ParallelFor(dat_crse.box(), [=] AMREX_GPU_DEVICE (int i, int j, int k)
466 {
467 int i_f = (i == -1) ? -1 : rx*i;
468 if (j == -1) {
469 fine_arr(i_f,j,k) = crse_arr(i,j,k);
470 } else {
471 fine_arr(i_f,ry*j ,k) = crse_arr(i,j,k);
472 if (j < bhi.y) {
473 for (int n = 1; n < ry; n++) {
474 fine_arr(i_f,ry*j+n,k) = crse_arr(i,j,k) + n * yfac * (crse_arr(i,j+1,k) - crse_arr(i,j,k));
475 }
476 }
477 }
478 });
479 } else {
480 amrex::ParallelFor(dat_crse.box(), [=] AMREX_GPU_DEVICE (int i, int j, int k)
481 {
482 int i_f = (i == -1) ? -1 : rx*i;
483 if (j == -1) {
484 fine_arr(i_f,j,k) = crse_arr(i,j,k);
485 } else {
486 fine_arr(i_f,ry*j,k) = crse_arr(i,j,k);
487 if (j < bhi.y) {
488 for (int n = 1; n < ry; n++) {
489 fine_arr(i_f,ry*j+n,k) = crse_arr(i,j,k) + n * yfac * (crse_arr(i,j+1,k) - crse_arr(i,j,k));
490 }
491 }
492 }
493 });
494 }
495 } else if (dat_crse.box().length(1) == 1) {
496 // Doing box on y-face
497 // amrex::Print() << "DOING INTERP ON YFACE " << dat_crse.box() << " " << dat_fine.box() << std::endl;
498 if (dat_crse.box().ixType()[0] == 0) {
499 amrex::ParallelFor(dat_crse.box(), [=] AMREX_GPU_DEVICE (int i, int j, int k)
500 {
501 int j_f = (j == -1) ? -1 : ry*j;
502 if (i == -1) {
503 fine_arr(i,j_f,k) = crse_arr(i,j,k);
504 } else {
505 fine_arr(rx*i ,j_f,k) = crse_arr(i,j,k);
506 if (i < bhi.x) {
507 for (int n = 1; n < rx; n++) {
508 fine_arr(rx*i+n,j_f,k) = crse_arr(i,j,k) + n * xfac * (crse_arr(i+1,j,k) - crse_arr(i,j,k));
509 }
510 }
511 }
512 });
513 } else {
514 amrex::ParallelFor(dat_crse.box(), [=] AMREX_GPU_DEVICE (int i, int j, int k)
515 {
516 int j_f = (j == -1) ? -1 : ry*j;
517 if (i == -1) {
518 fine_arr(i,j_f,k) = crse_arr(i,j,k);
519 } else {
520 fine_arr(rx*i,j_f,k) = crse_arr(i,j,k);
521 if (i < bhi.x) {
522 for (int n = 1; n < rx; n++) {
523 fine_arr(rx*i+n,j_f,k) = crse_arr(i,j,k) + n * xfac * (crse_arr(i+1,j,k) - crse_arr(i,j,k));
524 }
525 }
526 }
527 });
528 }
529 } else {
530 amrex::Abort(" What am I doing here??");
531 }
532}
533
534#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.
amrex::FArrayBox yhi_dat_before
FArrayBox to store data at time_before at yhi boundary.
amrex::FArrayBox yhi_dat_after
FArrayBox to store data at time_after at yhi boundary.
amrex::GpuArray< bool, AMREX_SPACEDIM *2 > var_need_data
Array over boundaries indicating whether we need physical data for this variable.
void interp_fab(amrex::FArrayBox &dat_crse, amrex::FArrayBox &dat_fine)
spatially interpolate boundary fabs from level 0 data (read in) to current level
amrex::FArrayBox ylo_dat_before
FArrayBox to store data at time_before at ylo boundary.
amrex::FArrayBox xlo_dat_after
FArrayBox to store data at time_after at xlo boundary.
amrex::Vector< std::string > nc_var_names
Variable names that will be read from file.
amrex::FArrayBox ylo_dat_after
FArrayBox to store data at time_after at ylo boundary.
amrex::Real time_before
Time in ocean_times immediately before the last time interpolated to.
void Initialize()
Read in time array from file and allocate data arrays.
amrex::FArrayBox yhi_dat_interp
FArrayBox to store data at inteprolated time at yhi boundary *‍/.
amrex::FArrayBox xlo_dat_interp
FArrayBox to store data at inteprolated time at xlo boundary *‍/.
amrex::Vector< amrex::Real > bry_times
Time points in netcdf file.
std::string time_name
Field name for time series in netcdf file.
amrex::FArrayBox xhi_dat_interp
FArrayBox to store data at inteprolated time at xhi boundary *‍/.
amrex::IntVect index_type
Index type for field to fill.
int i_time_before
Time index immediately before the last time interpolated to.
std::string field_name
Field name in netcdf file.
amrex::FArrayBox ylo_dat_interp
FArrayBox to store data at inteprolated time at ylo boundary *‍/.
amrex::Vector< amrex::Geometry > m_geom
Geometry at all levels.
amrex::FArrayBox xhi_dat_after
FArrayBox to store data at time_after at xhi boundary.
void read_in_at_time(amrex::FArrayBox &fab_xlo, amrex::FArrayBox &fab_xhi, amrex::FArrayBox &fab_ylo, amrex::FArrayBox &fab_yhi, int itime)
Read in data from file at time index itime and fill into mf.
amrex::Vector< std::string > file_names
File name to read from.
amrex::Vector< int > file_itime_offset
Offset to access a particular time within its file.
amrex::FArrayBox xhi_dat_before
FArrayBox to store data at time_before at xhi boundary.
bool is2d
Whether the field we're reading in is 2d.
void update_interpolated_to_time(amrex::Real time)
Calculate interpolated values at time, reading in data as necessary.
int m_lev
Level at which we are holding the boundary data.
int m_rx
Refinement ratios relative to level 0.
NCTimeSeriesBoundary(int a_lev, const amrex::Vector< amrex::Geometry > a_geom, const amrex::Vector< std::string > &a_file_name, const std::string a_field_name, const std::string a_time_name, const amrex::IntVect a_index_type, const amrex::GpuArray< bool, AMREX_SPACEDIM *2 > *a_var_need_data, bool a_is2d, int rx, int ry)
Constructor.
amrex::FArrayBox xlo_dat_before
FArrayBox to store data at time_before at xlo boundary.
amrex::Real time_after
Time in ocean_times immediately after the last time interpolated to.
amrex::Vector< int > file_for_time
File index to access a particular time.
NDArray is the datatype designed to hold any data, including scalars, multidimensional arrays,...