REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_NCTimeSeries.H
Go to the documentation of this file.
1#ifndef _REMORA_NCTIMESEIRES_H_
2#define _REMORA_NCTIMESEIRES_H_
3
4#ifdef REMORA_USE_NETCDF
5
6#include <string>
7#include <memory>
8
9#include <AMReX_AmrCore.H>
10
11/** \brief A class to hold and interpolate time series data read from a NetCDF file
12 *
13 * The class only ever holds the two time steps necessary to interpoplate to the current
14 * simulation time.
15 */
17{
18 public:
19 /** \brief Constructor */
20 NCTimeSeries (const amrex::Vector<std::string>& a_file_names,
21 const std::string a_field_name,
22 const std::string a_time_name,
23 const amrex::Box& a_domain,
24 amrex::MultiFab* a_mf_var,
25 bool a_is2d, bool a_save_interpolated);
26
27 /** \brief Read in time array from file and allocate data arrays */
28 void Initialize ();
29
30 /** \brief Calculate interpolated values at time and fill data for level lev */
31 void update_interpolated_to_time (amrex::Real time, int lev,
32 amrex::MultiFab* mf_lev,
33 const amrex::Vector<amrex::Geometry>& geom,
34 const amrex::Vector<amrex::IntVect>& ref_ratio);
35
36 /** \brief Access interpolated data saved for a specific level */
37 const amrex::MultiFab* get_interpolated_mf (int lev) const;
38 private:
39 /** \brief Read in data from file at time index itime and fill into mf */
40 void read_in_at_time (amrex::MultiFab* mf, int itime);
41
42 /** \brief Build cumulative refinement ratio from level 0 to lev */
43 amrex::IntVect cumulative_ref_ratio (int lev,
44 const amrex::Vector<amrex::IntVect>& ref_ratio) const;
45
46 /// File names to read from
47 amrex::Vector<std::string> file_names;
48 /// Field name in netcdf file
49 std::string field_name;
50 /// Field name for time series in netcdf file
51 std::string time_name;
52 /// Domain
53 amrex::Box domain;
54
55 /// Whether the field we're reading in is 2d
56 bool is2d;
57 /** Whether to save interpolated results in mf_interpolated. If false,
58 * data will be calculated to mf_var */
60
61 /// Time points in netcdf file
62 amrex::Vector<amrex::Real> ocean_times;
63 /// File index to access a particular time
64 amrex::Vector<int> file_for_time;
65 /// Offset to access a particular time within its file
66 amrex::Vector<int> file_itime_offset;
67 /// Time index immediately before the last time interpolated to
69 /// Time in ocean_times immediately before the last time interpolated to
70 amrex::Real time_before;
71 /// Time in ocean_times immediately after the last time interpolated to
72 amrex::Real time_after;
73
74 /// Multifab to store data at time_before
75 amrex::MultiFab* mf_before;
76 /// Multifab to store data at time_after
77 amrex::MultiFab* mf_after;
78 /// Multifab storing temporally interpolated data on level 0
79 amrex::MultiFab* mf_interp_lev0;
80 /// Interpolated data on each requested AMR level if save_interpolated=true
81 amrex::Vector<std::unique_ptr<amrex::MultiFab>> mf_interpolated_lev;
82 /** Pointer to REMORA data that corresponds to the variable being interpolated.
83 * Filled by update_interpolated_to_time if save_interpolated==false. Otherwise,
84 * just used for getting box array, nodality, and distribution mapping */
85 amrex::MultiFab* mf_var;
86};
87
88#endif // REMORA_USE_NETCDF
89#endif //_REMORA_NCTIMESEIRES_H_
A class to hold and interpolate time series data read from a NetCDF file.
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::Vector< int > file_itime_offset
Offset to access a particular time within its file.
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::Vector< std::string > file_names
File names to read from.
amrex::MultiFab * mf_after
Multifab to store data at time_after.
std::string time_name
Field name for time series in netcdf file.
amrex::Vector< int > file_for_time
File index to access a particular time.
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