22 #include "libmesh/equation_systems.h" 23 #include "libmesh/mesh_function.h" 24 #include "libmesh/numeric_vector.h" 25 #include "libmesh/nonlinear_implicit_system.h" 26 #include "libmesh/transient_system.h" 27 #include "libmesh/parallel_mesh.h" 28 #include "libmesh/serial_mesh.h" 29 #include "libmesh/exodusII_io.h" 30 #include "libmesh/exodusII_io_helper.h" 31 #include "libmesh/enum_xdr_mode.h" 32 #include "libmesh/string_to_enum.h" 42 "mesh",
"The name of the mesh file (must be xda/xdr or exodusII file).");
43 params.
addParam<std::vector<std::string>>(
45 std::vector<std::string>(),
46 "The name of the nodal and elemental variables from the file you want to use for values");
52 "The name of the file holding the equation system info in xda/xdr format (xda/xdr only).");
56 "The name of the system to pull values out of (xda/xdr only). The default name for the " 57 "nonlinear system is 'nl0', auxiliary system is 'aux0'");
60 params.
addParam<std::string>(
"timestep",
61 "Index of the single timestep used or \"LATEST\" for " 62 "the last timestep (exodusII only). If not supplied, " 63 "time interpolation will occur.");
68 "Specifies the order of the nodal solution data.");
72 "scale", std::vector<Real>(LIBMESH_DIM, 1),
"Scale factor for points in the simulation");
73 params.
addParam<std::vector<Real>>(
"scale_multiplier",
74 std::vector<Real>(LIBMESH_DIM, 1),
75 "Scale multiplying factor for points in the simulation");
76 params.
addParam<std::vector<Real>>(
"translation",
77 std::vector<Real>(LIBMESH_DIM, 0),
78 "Translation factors for x,y,z coordinates of the simulation");
81 "Vector about which to rotate points of the simulation.");
85 "Anticlockwise rotation angle (in degrees) to use for rotation about rotation0_vector.");
88 "Vector about which to rotate points of the simulation.");
92 "Anticlockwise rotation angle (in degrees) to use for rotation about rotation1_vector.");
96 "rotation0 translation scale rotation1 scale_multiplier",
"translation scale");
98 "transformation_order",
99 default_transformation_order,
100 "The order to perform the operations in. Define R0 to be the rotation matrix encoded by " 101 "rotation0_vector and rotation0_angle. Similarly for R1. Denote the scale by s, the " 102 "scale_multiplier by m, and the translation by t. Then, given a point x in the simulation, " 103 "if transformation_order = 'rotation0 scale_multiplier translation scale rotation1' then " 104 "form p = R1*(R0*x*m - t)/s. Then the values provided by the SolutionUserObjectBase at " 106 "in the simulation are the variable values at point p in the mesh.");
107 params.
addParamNamesToGroup(
"scale scale_multiplier translation rotation0_vector rotation0_angle " 108 "rotation1_angle transformation_order",
109 "Coordinate system transformation");
119 _file_type(
MooseEnum(
"xda=0 exodusII=1 xdr=2")),
120 _mesh_file(getParam<MeshFileName>(
"mesh")),
121 _es_file(getParam<FileName>(
"es")),
122 _system_name(getParam<
std::string>(
"system")),
123 _system_variables(getParam<
std::vector<
std::string>>(
"system_variables")),
124 _exodus_time_index(-1),
125 _interpolate_times(false),
128 _interpolation_time(0.0),
129 _interpolation_factor(0.0),
130 _exodus_times(nullptr),
133 _nodal_variable_order(getParam<
MooseEnum>(
"nodal_variable_order")),
134 _scale(getParam<
std::vector<
Real>>(
"scale")),
135 _scale_multiplier(getParam<
std::vector<
Real>>(
"scale_multiplier")),
136 _translation(getParam<
std::vector<
Real>>(
"translation")),
138 _rotation0_angle(getParam<
Real>(
"rotation0_angle")),
141 _rotation1_angle(getParam<
Real>(
"rotation1_angle")),
143 _transformation_order(getParam<
MultiMooseEnum>(
"transformation_order")),
171 if (
isParamValid(
"timestep") && getParam<std::string>(
"timestep") ==
"-1")
172 mooseError(
"A \"timestep\" of -1 is no longer supported for interpolation. Instead simply " 173 "remove this parameter altogether for interpolation");
178 "The mesh has second-order elements, be sure to set 'nodal_variable_order' if needed.");
185 paramError(
"es",
"Equation system file (.xda or .xdr) should have been specified");
195 _es = std::make_unique<EquationSystems>(*_mesh);
201 EquationSystems::READ_HEADER | EquationSystems::READ_DATA |
202 EquationSystems::READ_ADDITIONAL_DATA);
208 EquationSystems::READ_HEADER | EquationSystems::READ_DATA |
209 EquationSystems::READ_ADDITIONAL_DATA);
213 mooseError(
"Failed to determine proper read method for XDA/XDR equation system file: ",
229 _exodusII_io = std::make_unique<libMesh::ExodusII_IO>(*_mesh);
236 std::string s_timestep = getParam<std::string>(
"timestep");
238 if (s_timestep ==
"LATEST")
242 std::istringstream ss(s_timestep);
244 mooseError(
"Invalid value passed as \"timestep\". Expected \"LATEST\" or a valid integer " 257 if (num_exo_times == 0)
258 mooseError(
"In SolutionUserObjectBase, exodus file contains no timesteps.");
261 if (dynamic_cast<DistributedMesh *>(
_mesh.get()))
263 _mesh->allow_renumbering(
true);
264 _mesh->prepare_for_use();
268 _mesh->allow_renumbering(
false);
269 _mesh->prepare_for_use();
273 _es = std::make_unique<EquationSystems>(*_mesh);
278 const std::vector<std::string> & all_nodal(
_exodusII_io->get_nodal_var_names());
279 const std::vector<std::string> & all_elemental(
_exodusII_io->get_elem_var_names());
280 const std::vector<std::string> & all_scalar(
_exodusII_io->get_global_var_names());
293 if (
std::find(all_nodal.begin(), all_nodal.end(), var_name) != all_nodal.end())
295 else if (
std::find(all_elemental.begin(), all_elemental.end(), var_name) !=
298 else if (
std::find(all_scalar.begin(), all_scalar.end(), var_name) != all_scalar.end())
303 paramError(
"system_variables",
"Variable '" + var_name +
"' was not found in Exodus file");
311 for (
auto var_name : all_scalar)
338 _es2 = std::make_unique<EquationSystems>(*_mesh);
390 mooseError(
"In SolutionUserObjectBase, timestep = ",
392 ", but there are only ",
423 mooseAssert(sys_node.n_dofs(sys_num, var_num) > 0,
424 "Variable " << var_name <<
" has no DoFs on node " << sys_node.id());
441 mooseAssert(sys_elem.n_dofs(sys_num, var_num) > 0,
442 "Variable " << var_name <<
" has no DoFs on element " << sys_elem.id());
517 "In SolutionUserObjectBase, invalid file type: only .xda, .xdr, .exo and .e supported");
527 std::vector<unsigned int> var_nums;
533 for (
const auto & var_num : var_nums)
632 mooseError(
"In SolutionUserObjectBase, getTimeInterpolationData only applicable for exodusII " 650 for (
int i = 0; i < num_exo_times - 1; ++i)
657 (solution_time - (*_exodus_times)[i]) / ((*
_exodus_times)[i + 1] - (*_exodus_times)[i]);
660 else if (i == num_exo_times - 2)
670 bool indices_modified(
false);
673 indices_modified =
true;
675 return indices_modified;
684 mooseError(
"Value requested for nonexistent variable '",
688 "' SolutionUserObjectBase.\nSystem selected: ",
690 "\nAvailable variables:\n",
698 const std::string & var_name,
700 const std::set<subdomain_id_type> * subdomain_ids)
const 712 if (weighting_type == 1 ||
714 return pointValue(t, p, var_name, subdomain_ids);
718 switch (weighting_type)
723 for (
auto & v : values)
725 return average /
Real(values.size());
731 for (
auto & v : values)
732 if (v.first->id() < smallest_elem_id)
734 smallest_elem_id = v.first->id();
735 smallest_elem_id_value = v.second;
737 return smallest_elem_id_value;
741 Real largest_elem_id_value = std::numeric_limits<Real>::lowest();
743 for (
auto & v : values)
744 if (v.first->id() > largest_elem_id)
746 largest_elem_id = v.first->id();
747 largest_elem_id_value = v.second;
749 return largest_elem_id_value;
753 mooseError(
"SolutionUserObjectBase::pointValueWrapper reaches line that it should not be able to " 761 const std::string & var_name,
762 const std::set<subdomain_id_type> * subdomain_ids)
const 765 return pointValue(t, p, local_var_index, subdomain_ids);
771 const unsigned int local_var_index,
772 const std::set<subdomain_id_type> * subdomain_ids)
const 802 "Time passed into value() must match time at last call to timestepSetup()");
810 std::map<const Elem *, Real>
814 const std::string & var_name,
815 const std::set<subdomain_id_type> * subdomain_ids)
const 821 std::map<const Elem *, Real>
823 Real libmesh_dbg_var(t),
825 const unsigned int local_var_index,
826 const std::set<subdomain_id_type> * subdomain_ids)
const 847 std::map<const Elem *, Real> map =
854 "Time passed into value() must match time at last call to timestepSetup()");
857 if (map.size() != map2.size())
859 "In SolutionUserObjectBase::discontinuousPointValue map and map2 have different size");
864 if (map2.find(k.first) == map2.end())
866 "In SolutionUserObjectBase::discontinuousPointValue map and map2 have differing keys");
868 Real val2 = map2[k.first];
880 const std::string & var_name,
882 const std::set<subdomain_id_type> * subdomain_ids)
const 885 if (weighting_type == 1)
889 std::map<const Elem *, RealGradient> values =
891 switch (weighting_type)
896 for (
auto & v : values)
898 return average /
Real(values.size());
904 for (
auto & v : values)
905 if (v.first->id() < smallest_elem_id)
907 smallest_elem_id = v.first->id();
908 smallest_elem_id_value = v.second;
910 return smallest_elem_id_value;
916 for (
auto & v : values)
917 if (v.first->id() > largest_elem_id)
919 largest_elem_id = v.first->id();
920 largest_elem_id_value = v.second;
922 return largest_elem_id_value;
926 mooseError(
"SolutionUserObjectBase::pointValueGradientWrapper reaches line that it should not be " 934 const std::string & var_name,
935 const std::set<subdomain_id_type> * subdomain_ids)
const 944 const unsigned int local_var_index,
945 const std::set<subdomain_id_type> * subdomain_ids)
const 972 "Time passed into value() must match time at last call to timestepSetup()");
980 std::map<const Elem *, RealGradient>
984 const std::string & var_name,
985 const std::set<subdomain_id_type> * subdomain_ids)
const 991 std::map<const Elem *, RealGradient>
993 Real libmesh_dbg_var(t),
995 const unsigned int local_var_index,
996 const std::set<subdomain_id_type> * subdomain_ids)
const 1017 std::map<const Elem *, RealGradient> map =
1024 "Time passed into value() must match time at last call to timestepSetup()");
1025 std::map<const Elem *, RealGradient> map2 =
1028 if (map.size() != map2.size())
1030 "In SolutionUserObjectBase::discontinuousPointValue map and map2 have different size");
1033 for (
auto & k : map)
1035 if (map2.find(k.first) == map2.end())
1037 "In SolutionUserObjectBase::discontinuousPointValue map and map2 have differing keys");
1050 Real val = (*_serialized_solution)(dof_index);
1053 Real val2 = (*_serialized_solution2)(dof_index);
1061 const unsigned int local_var_index,
1062 unsigned int func_num,
1063 const std::set<subdomain_id_type> * subdomain_ids)
const 1079 (*_mesh_function)(p, 0.0, output, subdomain_ids);
1088 else if (func_num == 2)
1095 (*_mesh_function2)(p, 0.0, output, subdomain_ids);
1108 if (output.
size() == 0)
1110 std::ostringstream oss;
1112 mooseError(
"Failed to access the data for variable '",
1118 "' SolutionUserObjectBase");
1120 return output(local_var_index);
1123 std::map<const Elem *, Real>
1126 const unsigned int local_var_index,
1127 unsigned int func_num,
1128 const std::set<subdomain_id_type> * subdomain_ids)
const 1131 std::map<const Elem *, DenseVector<Number>> temporary_output;
1137 _mesh_function->discontinuous_value(p, 0.0, temporary_output, subdomain_ids);
1140 else if (func_num == 2)
1141 _mesh_function2->discontinuous_value(p, 0.0, temporary_output, subdomain_ids);
1149 if (temporary_output.size() == 0)
1151 std::ostringstream oss;
1153 mooseError(
"Failed to access the data for variable '",
1159 "' SolutionUserObjectBase");
1163 std::map<const Elem *, Real> output;
1164 for (
auto & k : temporary_output)
1167 k.second.size() > local_var_index,
1168 "In SolutionUserObjectBase::evalMultiValuedMeshFunction variable with local_var_index " 1169 << local_var_index <<
" does not exist");
1170 output[k.first] = k.second(local_var_index);
1179 const unsigned int local_var_index,
1180 unsigned int func_num,
1181 const std::set<subdomain_id_type> * subdomain_ids)
const 1184 std::vector<Gradient> output;
1193 else if (func_num == 2)
1202 if (output.size() == 0)
1204 std::ostringstream oss;
1206 mooseError(
"Failed to access the data for variable '",
1212 "' SolutionUserObjectBase");
1214 return output[local_var_index];
1217 std::map<const Elem *, RealGradient>
1220 const unsigned int local_var_index,
1221 unsigned int func_num,
1222 const std::set<subdomain_id_type> * subdomain_ids)
const 1225 std::map<const Elem *, std::vector<Gradient>> temporary_output;
1231 _mesh_function->discontinuous_gradient(p, 0.0, temporary_output, subdomain_ids);
1234 else if (func_num == 2)
1235 _mesh_function2->discontinuous_gradient(p, 0.0, temporary_output, subdomain_ids);
1243 if (temporary_output.size() == 0)
1245 std::ostringstream oss;
1247 mooseError(
"Failed to access the data for variable '",
1253 "' SolutionUserObjectBase");
1257 std::map<const Elem *, RealGradient> output;
1258 for (
auto & k : temporary_output)
1261 k.second.size() > local_var_index,
1262 "In SolutionUserObjectBase::evalMultiValuedMeshFunction variable with local_var_index " 1263 << local_var_index <<
" does not exist");
1264 output[k.first] = k.second[local_var_index];
1270 const std::vector<std::string> &
1288 std::vector<dof_id_type> dofs;
1297 #ifdef LIBMESH_HAVE_EXODUS_API 1302 for (
const auto & it : id_to_block)
MultiMooseEnum _transformation_order
transformations (rotations, translation, scales) are performed in this order
void mooseInfo(Args &&... args) const
std::vector< std::string > _system_variables
A list of variables to extract from the read system.
static InputParameters validParams()
SolutionUserObjectBase(const InputParameters ¶meters)
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
KOKKOS_INLINE_FUNCTION const T * find(const T &target, const T *const begin, const T *const end)
Find a value in an array.
std::unique_ptr< libMesh::MeshFunction > _mesh_function2
Pointer to second libMesh::MeshFuntion, used for interpolation.
const libMesh::FEType & feType() const
Get the type of finite element object.
std::unique_ptr< libMesh::EquationSystems > _es
Pointer to the libMesh::EquationSystems object.
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sin(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tan
RealVectorValue _rotation0_vector
vector about which to rotate
const std::vector< Real > * _exodus_times
The times available in the ExodusII file.
virtual void initialSetup() override
Initialize the System and Mesh objects for the solution being read.
virtual Real solutionSampleTime()=0
Get the time at which to sample the solution.
void paramError(const std::string ¶m, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
static InputParameters validParams()
void readXda()
Method for reading XDA mesh and equation systems file(s) This method is called by the constructor whe...
static Threads::spin_mutex _solution_user_object_mutex
std::map< const Elem *, libMesh::RealGradient > discontinuousPointValueGradient(Real t, const Point &p, const std::string &var_name, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Returns the gradient at a specific location and variable for cases where the gradient is multivalued ...
TypeTensor< T > transpose() const
std::unique_ptr< NumericVector< Number > > _serialized_solution2
Pointer to second serial solution, used for interpolation.
std::vector< Real > _scale
Scale parameter.
std::vector< Real > _translation
Translation.
Real _rotation1_angle
angle (in degrees) which to rotate through about vector _rotation1_vector
unsigned int size() const
Return the number of active items in the MultiMooseEnum.
virtual dof_id_type maxElemId() const
Real _rotation0_angle
angle (in degrees) which to rotate through about vector _rotation0_vector
static constexpr std::size_t dim
This is the dimension of all vector and tensor datastructures used in MOOSE.
MooseEnum getSolutionFileType() const
Get the type of file that was read.
void get_all_variable_numbers(std::vector< unsigned int > &all_variable_numbers) const
int _exodus_time_index
Current ExodusII time index.
std::unique_ptr< libMesh::ExodusII_IO > _exodusII_io
Pointer to the libMesh::ExodusII used to read the files.
Real _interpolation_time
Time in the current simulation at which the solution interpolation was last updated.
std::unique_ptr< libMesh::EquationSystems > _es2
Pointer to second libMesh::EquationSystems object, used for interpolation.
const Parallel::Communicator & _communicator
std::map< const Elem *, Real > evalMultiValuedMeshFunction(const Point &p, const unsigned int local_var_index, unsigned int func_num, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
A wrapper method for calling the various MeshFunctions that calls the mesh function functionality for...
libMesh::RealGradient evalMeshFunctionGradient(const Point &p, const unsigned int local_var_index, unsigned int func_num, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
A wrapper method interfacing with the libMesh mesh function for evaluating the gradient.
std::unique_ptr< NumericVector< Number > > _serialized_solution
Pointer to the serial solution vector.
Point _cached_p
Cached points.
bool _initialized
True if initial_setup has executed.
const MeshBase & get_mesh() const
virtual const MooseVariableFieldBase & getVariable(const THREAD_ID tid, const std::string &var_name, Moose::VarKindType expected_var_type=Moose::VarKindType::VAR_ANY, Moose::VarFieldType expected_var_field_type=Moose::VarFieldType::VAR_FIELD_ANY) const override
Returns the variable reference for requested variable which must be of the expected_var_type (Nonline...
dof_id_type n_dofs() const
std::map< SubdomainName, SubdomainID > _block_name_to_id
Map from block ID to block names. Read from the ExodusII file.
void SCALAR_dof_indices(std::vector< dof_id_type > &di, const unsigned int vn, const bool old_dofs=false) const
auto max(const L &left, const R &right)
unsigned int variable_number(std::string_view var) const
virtual void timestepSetup() override
When reading ExodusII files, this will update the interpolation times.
DenseVector< Number > _cached_values2
virtual void execute() override
Execute method.
DenseVector< Number > _cached_values
Cached values.
std::unique_ptr< libMesh::MeshBase > _mesh
Pointer the libMesh::mesh object.
libMesh::System * _system
Pointer libMesh::System class storing the read solution.
bool updateExodusBracketingTimeIndices()
Updates the time indices to interpolate between for ExodusII data.
unsigned int number() const
RealVectorValue _rotation1_vector
vector about which to rotate
const std::string & name() const
Get the name of the class.
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template cos(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(cos
std::string formatString(std::string message, const std::string &prefix)
Add new lines and prefixes to a string for pretty display in output NOTE: This makes a copy of the st...
std::string _mesh_file
The XDA or ExodusII file that is being read.
std::unique_ptr< NumericVector< Number > > solution
std::vector< std::string > _scalar_variables
Stores names of scalar variables.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
RealTensorValue _r1
Rotation matrix that performs the "_rotation1_angle about rotation1_vector".
unsigned int add_variable(std::string_view var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
const std::string & variable_name(const unsigned int i) const
void readExodusII()
Method for reading an ExodusII file, which is called when 'file_type = exodusII is set in the input f...
std::string stringify(const T &t)
conversion to string
bool hasSecondOrderElements()
check if the mesh has SECOND order elements
virtual void initialize() override
Called before execute() is ever called so that data can be cleared.
Real evalMeshFunction(const Point &p, const unsigned int local_var_index, unsigned int func_num, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
A wrapper method for calling the various MeshFunctions used for reading the data. ...
void updateExodusTimeInterpolation()
Updates the times for interpolating ExodusII data.
Real scalarValue(Real t, const std::string &var_name) const
Returns a value of a global variable.
const MooseEnum _nodal_variable_order
Nodal variable order, used when reading in solution data.
libMesh::System * _system2
Pointer to a second libMesh::System object, used for interpolation.
GenericRealTensorValue< is_ad > rotVecToZ(GenericRealVectorValue< is_ad > vec)
provides a rotation matrix that will rotate the vector vec to the z axis (the "2" direction) ...
bool hasExtension(const std::string &filename, std::string ext, bool strip_exodus_ext)
Real directValue(const Node *node, const std::string &var_name) const
Return a value directly from a Node.
void readBlockIdMapFromExodusII()
Read block ID map from the ExodusII file.
const std::vector< std::string > & variableNames() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
FEProblemBase & _fe_problem
Reference to the FEProblemBase for this user object.
bool isVariableNodal(const std::string &var_name) const
std::map< std::string, unsigned int > _local_variable_index
Stores the local index need by MeshFunction.
std::map< const Elem *, libMesh::RealGradient > evalMultiValuedMeshFunctionGradient(const Point &p, const unsigned int local_var_index, unsigned int func_num, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
A wrapper method interfacing with the libMesh mesh function that calls the gradient functionality for...
Real pointValueWrapper(Real t, const Point &p, const std::string &var_name, const MooseEnum &weighting_type=weightingType(), const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Returns a value at a specific location and variable checking for multiple values and weighting these ...
int _exodus_index1
Time index 1, used for interpolation.
virtual const Elem & elem_ref(const dof_id_type i) const
std::vector< std::string > _nodal_variables
Stores names of nodal variables.
std::map< int, std::string > id_to_block_names
std::string _es_file
The XDA file that contians the EquationSystems data (xda only)
const THREAD_ID _tid
Thread ID of this postprocessor.
IntRange< T > make_range(T beg, T end)
virtual MooseMesh & mesh() override
static std::unique_ptr< NumericVector< Number > > build(const Parallel::Communicator &comm, SolverPackage solver_package=libMesh::default_solver_package(), ParallelType parallel_type=AUTOMATIC)
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type and optionally a file path to the top-level block p...
RealTensorValue _r0
Rotation matrix that performs the "_rotation0_angle about rotation0_vector".
virtual unsigned int size() const override final
std::vector< std::string > _elemental_variables
Stores names of elemental variables.
virtual const Node & node_ref(const dof_id_type i) const
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type...
MooseEnum _file_type
File type to read (0 = xda; 1 = ExodusII)
libMesh::RealGradient pointValueGradientWrapper(Real t, const Point &p, const std::string &var_name, const MooseEnum &weighting_type=weightingType(), const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Returns the gradient at a specific location and variable checking for multiple values and weighting t...
std::vector< Real > _scale_multiplier
scale_multiplier parameter
std::unique_ptr< libMesh::MeshFunction > _mesh_function
Pointer the libMesh::MeshFunction object that the read data is stored.
bool _interpolate_times
Flag for triggering interpolation of ExodusII data.
bool checkFileReadable(const std::string &filename, bool check_line_endings, bool throw_on_unreadable, bool check_for_git_lfs_pointer)
std::string _system_name
The system name to extract from the XDA file (xda only)
const DofMap & get_dof_map() const
bool isParamSetByUser(const std::string &name) const
Test if the supplied parameter is set by a user, as opposed to not set or set to default.
std::map< const Elem *, Real > discontinuousPointValue(Real t, Point pt, const unsigned int local_var_index, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Returns a value at a specific location and variable for cases where the solution is multivalued at el...
Real pointValue(Real t, const Point &p, const unsigned int local_var_index, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Returns a value at a specific location and variable (see SolutionFunction)
unsigned int getLocalVarIndex(const std::string &var_name) const
Returns the local index for a given variable name.
int _exodus_index2
Time index 2, used for interpolation.
Real _interpolation_factor
Interpolation weight factor.
std::map< SubdomainID, SubdomainName > _block_id_to_name
Map from block names to block IDs. Read from the ExodusII file.
virtual void finalize() override
Finalize.
std::set< subdomain_id_type > _cached_subdomain_ids2
libMesh::RealGradient pointValueGradient(Real t, const Point &p, const std::string &var_name, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Returns the gradient at a specific location and variable (see SolutionFunction)
std::set< subdomain_id_type > _cached_subdomain_ids
Cached subdomain ids.