https://mooseframework.inl.gov
SolutionUserObjectBase.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 // MOOSE includes
13 #include "GeneralUserObject.h"
14 
15 // Forward declarations
16 namespace libMesh
17 {
18 class ExodusII_IO;
19 class EquationSystems;
20 class System;
21 class MeshFunction;
22 template <class T>
23 class NumericVector;
24 }
25 
31 {
32 public:
34 
36 
40  virtual Real solutionSampleTime() = 0;
41 
45  virtual void timestepSetup() override;
46 
52  unsigned int getLocalVarIndex(const std::string & var_name) const;
53 
66  const Point & p,
67  const std::string & var_name,
68  const MooseEnum & weighting_type = weightingType(),
69  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
70 
81  const Point & p,
82  const unsigned int local_var_index,
83  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
84 
95  const Point & p,
96  const std::string & var_name,
97  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
98 
111  std::map<const Elem *, Real>
113  Point pt,
114  const unsigned int local_var_index,
115  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
116 
129  std::map<const Elem *, Real>
131  const Point & p,
132  const std::string & var_name,
133  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
134 
148  const Point & p,
149  const std::string & var_name,
150  const MooseEnum & weighting_type = weightingType(),
151  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
152 
164  const Point & p,
165  const std::string & var_name,
166  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
167 
179  Point pt,
180  const unsigned int local_var_index,
181  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
182 
195  std::map<const Elem *, libMesh::RealGradient> discontinuousPointValueGradient(
196  Real t,
197  const Point & p,
198  const std::string & var_name,
199  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
200 
213  std::map<const Elem *, libMesh::RealGradient> discontinuousPointValueGradient(
214  Real t,
215  Point pt,
216  const unsigned int local_var_index,
217  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
218 
225  Real directValue(const Node * node, const std::string & var_name) const;
226 
233  Real directValue(const Elem * elem, const std::string & var_name) const;
234 
242  Real scalarValue(Real t, const std::string & var_name) const;
243 
244  // Required pure virtual function (not used)
245  virtual void initialize() override;
246 
247  // Required pure virtual function (not used)
248  virtual void finalize() override;
249 
250  // Required pure virtual function (not used)
251  virtual void execute() override;
252 
254  virtual void initialSetup() override;
255 
256  const std::vector<std::string> & variableNames() const;
257 
258  bool isVariableNodal(const std::string & var_name) const;
259 
261  {
262  return MooseEnum("found_first=1 average=2 smallest_element_id=4 largest_element_id=8",
263  "found_first");
264  }
265 
270  unsigned int getMeshFileDimension() const { return _mesh->spatial_dimension(); }
271 
275  const std::string getMeshFileName() const { return _mesh_file; }
276 
282  const std::map<SubdomainName, SubdomainID> & getBlockNamesToIds() const
283  {
284  return _block_name_to_id;
285  }
286 
292  const std::map<SubdomainID, SubdomainName> & getBlockIdsToNames() const
293  {
294  return _block_id_to_name;
295  }
296 
302 
303 protected:
309  void readXda();
310 
315  void readExodusII();
316 
324  virtual Real directValue(dof_id_type dof_index) const;
325 
330 
335 
343  Real evalMeshFunction(const Point & p,
344  const unsigned int local_var_index,
345  unsigned int func_num,
346  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
347 
356  std::map<const Elem *, Real>
357  evalMultiValuedMeshFunction(const Point & p,
358  const unsigned int local_var_index,
359  unsigned int func_num,
360  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
361 
370  evalMeshFunctionGradient(const Point & p,
371  const unsigned int local_var_index,
372  unsigned int func_num,
373  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
374 
384  std::map<const Elem *, libMesh::RealGradient> evalMultiValuedMeshFunctionGradient(
385  const Point & p,
386  const unsigned int local_var_index,
387  unsigned int func_num,
388  const std::set<subdomain_id_type> * subdomain_ids = nullptr) const;
389 
394 
397 
399  std::string _mesh_file;
400 
402  std::string _es_file;
403 
405  std::string _system_name;
406 
408  std::vector<std::string> _system_variables;
409 
411  std::map<std::string, unsigned int> _local_variable_index;
412 
414  std::vector<std::string> _nodal_variables;
415 
417  std::vector<std::string> _elemental_variables;
418 
420  std::vector<std::string> _scalar_variables;
421 
424 
427 
429  std::unique_ptr<libMesh::MeshBase> _mesh;
430 
432  std::unique_ptr<libMesh::EquationSystems> _es;
433 
436 
438  std::unique_ptr<libMesh::MeshFunction> _mesh_function;
439 
441  std::unique_ptr<libMesh::ExodusII_IO> _exodusII_io;
442 
444  std::unique_ptr<NumericVector<Number>> _serialized_solution;
445 
447  std::unique_ptr<libMesh::EquationSystems> _es2;
448 
451 
453  std::unique_ptr<libMesh::MeshFunction> _mesh_function2;
454 
456  std::unique_ptr<NumericVector<Number>> _serialized_solution2;
457 
460 
463 
465  const std::vector<Real> * _exodus_times;
466 
469 
472 
475 
477  std::vector<Real> _scale;
478 
480  std::vector<Real> _scale_multiplier;
481 
483  std::vector<Real> _translation;
484 
487 
490 
493 
496 
499 
502 
505 
508 
510  std::map<SubdomainName, SubdomainID> _block_name_to_id;
511 
513  std::map<SubdomainID, SubdomainName> _block_id_to_name;
514 
515 private:
517 };
MultiMooseEnum _transformation_order
transformations (rotations, translation, scales) are performed in this order
std::vector< std::string > _system_variables
A list of variables to extract from the read system.
static InputParameters validParams()
SolutionUserObjectBase(const InputParameters &parameters)
std::unique_ptr< libMesh::MeshFunction > _mesh_function2
Pointer to second libMesh::MeshFuntion, used for interpolation.
std::unique_ptr< libMesh::EquationSystems > _es
Pointer to the libMesh::EquationSystems object.
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 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 ...
User object that reads an existing solution from an input file and uses it in the current simulation...
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
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
Real _rotation0_angle
angle (in degrees) which to rotate through about vector _rotation0_vector
static MooseEnum weightingType()
MooseEnum getSolutionFileType() const
Get the type of file that was read.
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.
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...
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
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.
bool _initialized
True if initial_setup has executed.
std::map< SubdomainName, SubdomainID > _block_name_to_id
Map from block ID to block names. Read from the ExodusII file.
virtual void timestepSetup() override
When reading ExodusII files, this will update the interpolation times.
virtual void execute() override
Execute method.
const std::string getMeshFileName() const
Return the name of the mesh file this object read the solution from.
std::unique_ptr< libMesh::MeshBase > _mesh
Pointer the libMesh::mesh object.
TensorValue< Real > RealTensorValue
libMesh::System * _system
Pointer libMesh::System class storing the read solution.
bool updateExodusBracketingTimeIndices()
Updates the time indices to interpolate between for ExodusII data.
RealVectorValue _rotation1_vector
vector about which to rotate
std::string _mesh_file
The XDA or ExodusII file that is being read.
std::vector< std::string > _scalar_variables
Stores names of scalar variables.
unsigned int getMeshFileDimension() const
Return the spatial dimension of the mesh file.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
RealTensorValue _r1
Rotation matrix that performs the "_rotation1_angle about rotation1_vector".
void readExodusII()
Method for reading an ExodusII file, which is called when &#39;file_type = exodusII is set in the input f...
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.
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
const std::map< SubdomainName, SubdomainID > & getBlockNamesToIds() const
Get the map from block name to block ID.
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.
std::vector< std::string > _nodal_variables
Stores names of nodal variables.
std::string _es_file
The XDA file that contians the EquationSystems data (xda only)
RealTensorValue _r0
Rotation matrix that performs the "_rotation0_angle about rotation0_vector".
std::vector< std::string > _elemental_variables
Stores names of elemental variables.
const InputParameters & parameters() const
Get the parameters of the object.
const std::map< SubdomainID, SubdomainName > & getBlockIdsToNames() const
Get the map from block id to block name.
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.
std::string _system_name
The system name to extract from the XDA file (xda only)
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.
uint8_t dof_id_type
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)