PRISMS-PF Manual v3.0-pre
Loading...
Searching...
No Matches
variableContainer.h
1// This class permits the access of a subset of indexed fields and gives an
2// error if any non-allowed fields are requested
3#ifndef VARIBLECONTAINER_H
4#define VARIBLECONTAINER_H
5
6#include <deal.II/lac/vector.h>
7#include <deal.II/matrix_free/evaluation_flags.h>
8#include <deal.II/matrix_free/fe_evaluation.h>
9#include <deal.II/matrix_free/matrix_free.h>
10
11#include <boost/unordered_map.hpp>
12
13#include <core/model_variables.h>
14
15template <int dim, int degree, typename T>
17{
18public:
19#include <core/typeDefs.h>
20 // Constructors
21
22 // Standard contructor, used for most situations
23 variableContainer(const dealii::MatrixFree<dim, double> &data,
24 const std::vector<variable_info> &_varInfoList,
25 const std::vector<variable_info> &_varChangeInfoList);
26
27 variableContainer(const dealii::MatrixFree<dim, double> &data,
28 const std::vector<variable_info> &_varInfoList);
29 // Nonstandard constructor, used when only one index of "data" should be used,
30 // use with care!
31 variableContainer(const dealii::MatrixFree<dim, double> &data,
32 const std::vector<variable_info> &_varInfoList,
33 const unsigned int &fixed_index);
34
35 // Methods to get the value/grad/hess in the residual method (this is how the
36 // user gets these values in equations.h)
37 T
38 get_scalar_value(unsigned int global_variable_index) const;
39 dealii::Tensor<1, dim, T>
40 get_scalar_gradient(unsigned int global_variable_index) const;
41 dealii::Tensor<2, dim, T>
42 get_scalar_hessian(unsigned int global_variable_index) const;
43 dealii::Tensor<1, dim, T>
44 get_vector_value(unsigned int global_variable_index) const;
45 dealii::Tensor<2, dim, T>
46 get_vector_gradient(unsigned int global_variable_index) const;
47 dealii::Tensor<3, dim, T>
48 get_vector_hessian(unsigned int global_variable_index) const;
49
50 T
51 get_change_in_scalar_value(unsigned int global_variable_index) const;
52 dealii::Tensor<1, dim, T>
53 get_change_in_scalar_gradient(unsigned int global_variable_index) const;
54 dealii::Tensor<2, dim, T>
55 get_change_in_scalar_hessian(unsigned int global_variable_index) const;
56 dealii::Tensor<1, dim, T>
57 get_change_in_vector_value(unsigned int global_variable_index) const;
58 dealii::Tensor<2, dim, T>
59 get_change_in_vector_gradient(unsigned int global_variable_index) const;
60 dealii::Tensor<3, dim, T>
61 get_change_in_vector_hessian(unsigned int global_variable_index) const;
62
63 // Methods to set the value residual and the gradient residual (this is how
64 // the user sets these values in equations.h)
65 void
66 set_scalar_value_term_RHS(unsigned int global_variable_index, T val);
67 void
68 set_scalar_gradient_term_RHS(unsigned int global_variable_index,
69 dealii::Tensor<1, dim, T> grad);
70 void
71 set_vector_value_term_RHS(unsigned int global_variable_index,
72 dealii::Tensor<1, dim, T> val);
73 void
74 set_vector_gradient_term_RHS(unsigned int global_variable_index,
75 dealii::Tensor<2, dim, T> grad);
76
77 void
78 set_scalar_value_term_LHS(unsigned int global_variable_index, T val);
79 void
80 set_scalar_gradient_term_LHS(unsigned int global_variable_index,
81 dealii::Tensor<1, dim, T> grad);
82 void
83 set_vector_value_term_LHS(unsigned int global_variable_index,
84 dealii::Tensor<1, dim, T> val);
85 void
86 set_vector_gradient_term_LHS(unsigned int global_variable_index,
87 dealii::Tensor<2, dim, T> grad);
88
89 // Initialize, read DOFs, and set evaulation flags for each variable
90 void
91 reinit_and_eval(
92 const std::vector<dealii::LinearAlgebra::distributed::Vector<double> *> &src,
93 unsigned int cell);
94 void
96 const dealii::LinearAlgebra::distributed::Vector<double> &src,
97 unsigned int cell,
98 unsigned int var_being_solved);
99
100 // Only initialize the FEEvaluation object for each variable (used for
101 // post-processing)
102 void
103 reinit(unsigned int cell);
104
105 // Integrate the residuals and distribute from local to global
106 void
107 integrate_and_distribute(
108 std::vector<dealii::LinearAlgebra::distributed::Vector<double> *> &dst);
109 void
110 integrate_and_distribute_change_in_solution_LHS(
111 dealii::LinearAlgebra::distributed::Vector<double> &dst,
112 const unsigned int var_being_solved);
113
114 // The quadrature point index, a method to get the number of quadrature points
115 // per cell, and a method to get the xyz coordinates for the quadrature point
116 unsigned int q_point;
117
118 unsigned int
119 get_num_q_points() const;
120
121 dealii::Point<dim, T>
122 get_q_point_location() const;
123
124private:
125 // Vectors of the actual FEEvaluation objects for each active variable, split
126 // into scalar variables and vector variables for type reasons
127 using scalar_FEEval = dealii::FEEvaluation<dim, degree, degree + 1, 1, double>;
128 using vector_FEEval = dealii::FEEvaluation<dim, degree, degree + 1, dim, double>;
129
130 boost::unordered_map<unsigned int, std::unique_ptr<scalar_FEEval>> scalar_vars_map;
131 boost::unordered_map<unsigned int, std::unique_ptr<vector_FEEval>> vector_vars_map;
132
133 boost::unordered_map<unsigned int, std::unique_ptr<scalar_FEEval>>
134 scalar_change_in_vars_map;
135 boost::unordered_map<unsigned int, std::unique_ptr<vector_FEEval>>
136 vector_change_in_vars_map;
137
138 // Object containing some information about each variable (indices, whether
139 // the val/grad/hess is needed, etc)
140 std::vector<variable_info> varInfoList;
141 std::vector<variable_info> varChangeInfoList;
142
143 // The number of variables
144 unsigned int num_var;
145};
146
147#endif
Definition variableContainer.h:17
void reinit_and_eval_change_in_solution(const dealii::LinearAlgebra::distributed::Vector< double > &src, unsigned int cell, unsigned int var_being_solved)
Definition variableContainer.cc:213