PRISMS-PF Manual v3.0-pre
All Classes Functions Variables Enumerations Pages
variable_container.h
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#ifndef variable_container_h
5#define variable_container_h
6
7#include <deal.II/lac/vector.h>
8#include <deal.II/matrix_free/evaluation_flags.h>
9#include <deal.II/matrix_free/fe_evaluation.h>
10#include <deal.II/matrix_free/matrix_free.h>
11
12#include <prismspf/config.h>
13#include <prismspf/core/exceptions.h>
14#include <prismspf/core/type_enums.h>
15#include <prismspf/core/variable_attributes.h>
16
17PRISMS_PF_BEGIN_NAMESPACE
18
28template <int dim, int degree, typename number>
30{
31public:
32 using VectorType = dealii::LinearAlgebra::distributed::Vector<number>;
33 using value_type = number;
34 using size_type = dealii::VectorizedArray<number>;
35
39 variableContainer(const dealii::MatrixFree<dim, number> &data,
40 const std::map<unsigned int, variableAttributes> &_subset_attributes,
41 const std::unordered_map<std::pair<unsigned int, dependencyType>,
42 unsigned int,
43 pairHash> &_global_to_local_solution,
44 const solveType &_solve_type);
45
49 size_type
50 get_scalar_value(unsigned int global_variable_index,
51 dependencyType dependency_type = dependencyType::NORMAL) const;
52
56 dealii::Tensor<1, dim, size_type>
57 get_scalar_gradient(unsigned int global_variable_index,
58 dependencyType dependency_type = dependencyType::NORMAL) const;
59
63 dealii::Tensor<2, dim, size_type>
64 get_scalar_hessian(unsigned int global_variable_index,
65 dependencyType dependency_type = dependencyType::NORMAL) const;
66
70 dealii::Tensor<1, dim, size_type>
72 unsigned int global_variable_index,
73 dependencyType dependency_type = dependencyType::NORMAL) const;
74
78 size_type
79 get_scalar_laplacian(unsigned int global_variable_index,
80 dependencyType dependency_type = dependencyType::NORMAL) const;
81
85 dealii::Tensor<1, dim, size_type>
86 get_vector_value(unsigned int global_variable_index,
87 dependencyType dependency_type = dependencyType::NORMAL) const;
88
92 dealii::Tensor<2, dim, size_type>
93 get_vector_gradient(unsigned int global_variable_index,
94 dependencyType dependency_type = dependencyType::NORMAL) const;
95
99 dealii::Tensor<3, dim, size_type>
100 get_vector_hessian(unsigned int global_variable_index,
101 dependencyType dependency_type = dependencyType::NORMAL) const;
102
106 dealii::Tensor<2, dim, size_type>
108 unsigned int global_variable_index,
109 dependencyType dependency_type = dependencyType::NORMAL) const;
110
114 dealii::Tensor<1, dim, size_type>
115 get_vector_laplacian(unsigned int global_variable_index,
116 dependencyType dependency_type = dependencyType::NORMAL) const;
117
121 size_type
122 get_vector_divergence(unsigned int global_variable_index,
123 dependencyType dependency_type = dependencyType::NORMAL) const;
124
128 dealii::Tensor<2, dim, size_type>
130 unsigned int global_variable_index,
131 dependencyType dependency_type = dependencyType::NORMAL) const;
132
138 dealii::Tensor<1, (dim == 2 ? 1 : dim), size_type>
139 get_vector_curl(unsigned int global_variable_index,
140 dependencyType dependency_type = dependencyType::NORMAL) const;
141
145 void
146 set_scalar_value_term(const unsigned int &global_variable_index,
147 const size_type &val,
148 const dependencyType &dependency_type = dependencyType::NORMAL);
149
153 void
155 const unsigned int &global_variable_index,
156 const dealii::Tensor<1, dim, size_type> &grad,
157 const dependencyType &dependency_type = dependencyType::NORMAL);
158
162 void
163 set_vector_value_term(const unsigned int &global_variable_index,
164 const dealii::Tensor<1, dim, size_type> &val,
165 const dependencyType &dependency_type = dependencyType::NORMAL);
166
170 void
172 const unsigned int &global_variable_index,
173 const dealii::Tensor<2, dim, size_type> &grad,
174 const dependencyType &dependency_type = dependencyType::NORMAL);
175
180 void
182 const std::function<void(variableContainer &, const dealii::Point<dim, size_type> &)>
183 &func,
184 std::vector<VectorType *> &dst,
185 const std::vector<VectorType *> &src,
186 const std::pair<unsigned int, unsigned int> &cell_range);
187
192 void
194 const std::function<void(variableContainer &, const dealii::Point<dim, size_type> &)>
195 &func,
196 VectorType &dst,
197 const std::vector<VectorType *> &src,
198 const std::pair<unsigned int, unsigned int> &cell_range);
199
204 void
206 const std::function<void(variableContainer &, const dealii::Point<dim, size_type> &)>
207 &func,
208 VectorType &dst,
209 const VectorType &src,
210 const std::vector<VectorType *> &src_subset,
211 const std::pair<unsigned int, unsigned int> &cell_range);
212
216 void
218 const std::function<void(variableContainer &, const dealii::Point<dim, size_type> &)>
219 &func,
220 VectorType &dst,
221 const std::vector<VectorType *> &src_subset,
222 const std::pair<unsigned int, unsigned int> &cell_range);
223
224private:
225 using scalar_FEEval = dealii::FEEvaluation<dim, degree, degree + 1, 1, number>;
226 using vector_FEEval = dealii::FEEvaluation<dim, degree, degree + 1, dim, number>;
227
231 void
232 scalar_FEEval_exists(const unsigned int &dependency_index,
233 const dependencyType &dependency_type) const;
234
238 void
239 vector_FEEval_exists(const unsigned int &dependency_index,
240 const dependencyType &dependency_type) const;
241
246 void
247 access_valid(const unsigned int &dependency_index,
248 const dependencyType &dependency_type,
249 const dealii::EvaluationFlags::EvaluationFlags &flag) const;
250
254 void
255 submission_valid(const dependencyType &dependency_type) const;
256
260 [[nodiscard]] unsigned int
261 get_n_q_points() const;
262
266 [[nodiscard]] dealii::Point<dim, size_type>
267 get_q_point_location() const;
268
272 void
273 reinit_and_eval(const std::vector<VectorType *> &src, unsigned int cell);
274
278 void
279 reinit_and_eval(const VectorType &src, unsigned int cell);
280
284 void
285 reinit(unsigned int cell, const unsigned int &global_variable_index);
286
290 void
291 read_dof_values(const std::vector<VectorType *> &src, unsigned int cell);
292
297 void
298 eval(const unsigned int &global_variable_index);
299
303 void
304 integrate_and_distribute(std::vector<VectorType *> &dst);
305
309 void
310 integrate_and_distribute(VectorType &dst);
311
315 void
316 integrate(const unsigned int &global_variable_index);
317
322 std::map<unsigned int, std::map<dependencyType, std::unique_ptr<scalar_FEEval>>>
323 scalar_vars_map;
324
329 std::map<unsigned int, std::map<dependencyType, std::unique_ptr<vector_FEEval>>>
330 vector_vars_map;
331
335 const std::map<unsigned int, variableAttributes> &subset_attributes;
336
340 const std::unordered_map<std::pair<unsigned int, dependencyType>,
341 unsigned int,
342 pairHash> &global_to_local_solution;
343
349 std::unordered_map<std::pair<unsigned int, dependencyType>,
350 dealii::EvaluationFlags::EvaluationFlags,
351 pairHash>
352 src_eval_flags;
353
359 dealii::EvaluationFlags::EvaluationFlags dst_eval_flags =
360 dealii::EvaluationFlags::EvaluationFlags::nothing;
361
365 const solveType solve_type;
366
370 unsigned int q_point = 0;
371
375 unsigned int n_dofs_per_cell;
376
380 std::unique_ptr<dealii::AlignedVector<size_type>> scalar_diagonal;
381
385 std::unique_ptr<dealii::AlignedVector<dealii::Tensor<1, dim, size_type>>>
386 vector_diagonal;
387};
388
389PRISMS_PF_END_NAMESPACE
390
391#endif
This class permits the access of a subset of indexed fields and gives an error if any non-allowed fie...
Definition variable_container.h:30
dealii::Tensor< 3, dim, size_type > get_vector_hessian(unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const
Return the hessian of the specified vector field.
Definition variable_container.cc:1166
void eval_local_operator(const std::function< void(variableContainer &, const dealii::Point< dim, size_type > &)> &func, std::vector< VectorType * > &dst, const std::vector< VectorType * > &src, const std::pair< unsigned int, unsigned int > &cell_range)
Apply some operator function for a given cell range and source vector to some destination vector.
Definition variable_container.cc:78
size_type get_scalar_value(unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const
Return the value of the specified scalar field.
Definition variable_container.cc:1016
dealii::Tensor< 1, dim, size_type > get_vector_laplacian(unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const
Return the laplacian of the specified vector field.
Definition variable_container.cc:1227
void eval_local_diagonal(const std::function< void(variableContainer &, const dealii::Point< dim, size_type > &)> &func, VectorType &dst, const std::vector< VectorType * > &src_subset, const std::pair< unsigned int, unsigned int > &cell_range)
TODO (landinjm): Add comments.
Definition variable_container.cc:173
size_type get_vector_divergence(unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const
Return the divergence of the specified vector field.
Definition variable_container.cc:1257
dealii::Tensor< 1, dim, size_type > get_scalar_hessian_diagonal(unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const
Return the diagonal of the hessian of the specified scalar field.
Definition variable_container.cc:1070
void set_scalar_value_term(const unsigned int &global_variable_index, const size_type &val, const dependencyType &dependency_type=dependencyType::NORMAL)
Set the residual value of the specified scalar field.
Definition variable_container.cc:1322
dealii::Tensor< 1, dim, size_type > get_vector_value(unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const
Return the value of the specified vector field.
Definition variable_container.cc:1106
dealii::Tensor< 2, dim, size_type > get_vector_hessian_diagonal(unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const
Return the diagonal of the hessian of the specified vector field.
Definition variable_container.cc:1196
void set_vector_gradient_term(const unsigned int &global_variable_index, const dealii::Tensor< 2, dim, size_type > &grad, const dependencyType &dependency_type=dependencyType::NORMAL)
Set the residual gradient of the specified vector field.
Definition variable_container.cc:1374
size_type get_scalar_laplacian(unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const
Return the laplacian of the specified scalar field.
Definition variable_container.cc:1088
dealii::Tensor< 1, dim, size_type > get_scalar_gradient(unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const
Return the gradient of the specified scalar field.
Definition variable_container.cc:1034
dealii::Tensor< 2, dim, size_type > get_vector_gradient(unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const
Return the gradient of the specified vector field.
Definition variable_container.cc:1136
dealii::Tensor< 1,(dim==2 ? 1 :dim), size_type > get_vector_curl(unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const
Return the curl of the specified vector field. Note that this is dealii::VectorizedArray<number> type...
Definition variable_container.cc:1295
dealii::Tensor< 2, dim, size_type > get_scalar_hessian(unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const
Return the hessian of the specified scalar field.
Definition variable_container.cc:1052
void set_scalar_gradient_term(const unsigned int &global_variable_index, const dealii::Tensor< 1, dim, size_type > &grad, const dependencyType &dependency_type=dependencyType::NORMAL)
Set the residual gradient of the specified scalar field.
Definition variable_container.cc:1339
void set_vector_value_term(const unsigned int &global_variable_index, const dealii::Tensor< 1, dim, size_type > &val, const dependencyType &dependency_type=dependencyType::NORMAL)
Set the residual value of the specified vector field.
Definition variable_container.cc:1356
dealii::Tensor< 2, dim, size_type > get_vector_symmetric_gradient(unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const
Return the symmetric gradient of the specified vector field.
Definition variable_container.cc:1275
Simple hash function for pairs.
Definition variable_attributes.h:24