PRISMS-PF Manual
Loading...
Searching...
No Matches
invm_manager.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#pragma once
5
6#include <deal.II/base/vectorization.h>
7#include <deal.II/dofs/dof_handler.h>
8#include <deal.II/fe/fe_system.h>
9#include <deal.II/lac/affine_constraints.h>
10#include <deal.II/lac/la_parallel_vector.h>
11#include <deal.II/matrix_free/fe_evaluation.h>
12#include <deal.II/matrix_free/matrix_free.h>
13
18
19#include <prismspf/config.h>
20
22
26template <unsigned int dim, unsigned int degree, typename number>
28{
29public:
30 using ScalarValue = dealii::VectorizedArray<number>;
31 using VectorValue = dealii::Tensor<1, dim, ScalarValue>;
32
39 explicit InvMManager(const DoFManager<dim, degree> &dof_manager,
40 const ConstraintManager<dim, degree, number> &constraint_manager,
43 : num_levels(dof_manager.get_dof_handlers().size())
46 {
47 data.resize(num_levels);
48 jxw_scalar.resize(num_levels);
49 invm_scalar.resize(num_levels);
51 jxw_vector.resize(num_levels);
52 invm_vector.resize(num_levels);
54 reinit(dof_manager, constraint_manager);
55 }
56
57 void
58 reinit(const DoFManager<dim, degree> &dof_manager,
59 const ConstraintManager<dim, degree, number> &constraint_manager)
60 {
61 const std::vector<std::array<dealii::DoFHandler<dim>, 2>> &dof_handlers =
62 dof_manager.get_dof_handlers();
63 const std::vector<std::array<dealii::AffineConstraints<number>, 2>> &constraints =
64 constraint_manager.get_generic_constraints();
65 for (unsigned int i = 0; i < data.size(); ++i)
66 {
67 for (unsigned int rank = 0; rank < 2; ++rank)
68 {
70 dof_handlers[i][rank],
71 constraints[i][rank],
73 }
74 }
75 }
76
80 void
82 {
83 initialize();
85 {
87 }
89 {
91 }
92 }
93
101 get_invm(TensorRank rank, unsigned int relative_level) const
102 {
105 dealii::ExcInternalError("Requested invm that was not calculated"));
107 {
108 return invm_scalar[relative_level];
109 }
110 // else
111 {
112 return invm_vector[relative_level];
113 }
114 }
115
123 get_jxw(TensorRank rank, unsigned int relative_level) const
124 {
127 dealii::ExcInternalError("Requested invm that was not calculated"));
129 {
130 return jxw_scalar[relative_level];
131 }
132 // else
133 {
134 return jxw_vector[relative_level];
135 }
136 }
137
145 get_invm_sqrt(TensorRank rank, unsigned int relative_level) const
146 {
149 dealii::ExcInternalError("Requested invm that was not calculated"));
151 {
152 return invm_sqrt_scalar[relative_level];
153 }
154 // else
155 {
156 return invm_sqrt_vector[relative_level];
157 }
158 }
159
160 std::vector<const SolutionVector<number> *>
161 get_invm(const std::vector<FieldAttributes> &field_container,
162 const std::set<unsigned int> &field_indices,
163 unsigned int relative_level) const
164 {
165 const unsigned int num_blocks = field_indices.size();
166 std::vector<const SolutionVector<number> *> out;
167 out.reserve(num_blocks);
168 for (unsigned int field_index : field_indices)
169 {
170 out.push_back(&get_invm(field_container[field_index].field_type, relative_level));
171 }
172 return out;
173 }
174
175 std::vector<const SolutionVector<number> *>
176 get_jxw(const std::vector<FieldAttributes> &field_container,
177 const std::set<unsigned int> &field_indices,
178 unsigned int relative_level) const
179 {
180 const unsigned int num_blocks = field_indices.size();
181 std::vector<const SolutionVector<number> *> out;
182 out.reserve(num_blocks);
183 for (unsigned int field_index : field_indices)
184 {
185 out.push_back(&get_jxw(field_container[field_index].field_type, relative_level));
186 }
187 return out;
188 }
189
190 std::vector<const SolutionVector<number> *>
191 get_invm_sqrt(const std::vector<FieldAttributes> &field_container,
192 const std::set<unsigned int> &field_indices,
193 unsigned int relative_level) const
194 {
195 const unsigned int num_blocks = field_indices.size();
196 std::vector<const SolutionVector<number> *> out;
197 out.reserve(num_blocks);
198 for (unsigned int field_index : field_indices)
199 {
200 out.push_back(
201 &get_invm_sqrt(field_container[field_index].field_type, relative_level));
202 }
203 return out;
204 }
205
206private:
210 void
212 {
213 for (unsigned int i = 0; i < num_levels; ++i)
214 {
216 {
220 }
222 {
226 }
227 }
228 }
229
233 void
235 {
236 for (unsigned int i = 0; i < data.size(); ++i)
237 {
238 data[i][0].cell_loop(&InvMManager::compute_local_scalar, this, jxw_scalar[i], 0);
240 //
242 jxw_scalar[i].update_ghost_values();
243 invm_scalar[i].update_ghost_values();
244 invm_sqrt_scalar[i].update_ghost_values();
245 }
246 }
247
248 void
250 {
251 for (unsigned int i = 0; i < data.size(); ++i)
252 {
253 data[i][1].cell_loop(&InvMManager::compute_local_vector, this, jxw_vector[i], 0);
255 //
257 jxw_vector[i].update_ghost_values();
258 invm_vector[i].update_ghost_values();
259 invm_sqrt_vector[i].update_ghost_values();
260 }
261 }
262
263 void
264 compute_local_scalar(const MatrixFree<dim, number> &_data,
266 [[maybe_unused]] const int &src,
267 const std::pair<unsigned int, unsigned int> &cell_range) const
268 {
269 dealii::FEEvaluation<dim, degree, degree + 1, 1, number, ScalarValue> fe_eval(_data);
270 for (unsigned int cell = cell_range.first; cell < cell_range.second; ++cell)
271 {
272 fe_eval.reinit(cell);
273 for (unsigned int quad = 0; quad < SystemWide<dim, degree>::quadrature.size();
274 ++quad)
275 {
276 fe_eval.submit_value(dealii::make_vectorized_array<number>(1.0), quad);
277 }
278 fe_eval.integrate_scatter(dealii::EvaluationFlags::values, dst);
279 }
280 }
281
282 void
283 compute_local_vector(const MatrixFree<dim, number> &_data,
285 [[maybe_unused]] const int &src,
286 const std::pair<unsigned int, unsigned int> &cell_range) const
287 {
288 dealii::FEEvaluation<dim, degree, degree + 1, dim, number> fe_eval(_data);
289 for (unsigned int cell = cell_range.first; cell < cell_range.second; ++cell)
290 {
291 fe_eval.reinit(cell);
292 for (unsigned int quad = 0; quad < SystemWide<dim, degree>::quadrature.size();
293 ++quad)
294 {
295 fe_eval.submit_value(one, quad);
296 }
297 fe_eval.integrate_scatter(dealii::EvaluationFlags::values, dst);
298 }
299 }
300
301 void
303 {
304 src.update_ghost_values();
305 dst.update_ghost_values();
306 Assert(dst.size() == src.size(), dealii::ExcInternalError());
307 for (unsigned int i = 0; i < src.locally_owned_size(); ++i)
308 {
309 number src_el = src.local_element(i);
310 // for some reason, some elements are zero
311 dst.local_element(i) = src_el == 0 ? 0.0 : 1.0 / src_el;
312 }
313 }
314
315 void
317 {
318 using std::sqrt;
319 src.update_ghost_values();
320 dst.update_ghost_values();
321 Assert(dst.size() == src.size(), dealii::ExcInternalError());
322 for (unsigned int i = 0; i < src.locally_owned_size(); ++i)
323 {
324 dst.local_element(i) = sqrt(src.local_element(i));
325 }
326 }
327
331 std::vector<std::array<MatrixFree<dim, number>, 2>> data;
332
333 unsigned int num_levels;
334
335 bool calculate_scalar = false;
336 bool calculate_vector = false;
337
341 std::vector<SolutionVector<number>> jxw_scalar;
342 std::vector<SolutionVector<number>> jxw_vector;
343 std::vector<SolutionVector<number>> invm_scalar;
344 std::vector<SolutionVector<number>> invm_vector;
345 std::vector<SolutionVector<number>> invm_sqrt_scalar;
346 std::vector<SolutionVector<number>> invm_sqrt_vector;
347
348 inline static const VectorValue one = []()
349 {
351 for (unsigned int i = 0; i < dim; ++i)
352 {
353 one1[i] = 1.0;
354 }
355 return one1;
356 }();
357};
358
The class handles the generation and application of boundary conditions based on the user-inputs.
Definition constraint_manager.h:50
const std::vector< std::array< dealii::AffineConstraints< number >, 2 > > & get_generic_constraints() const
Getter function for the constraints.
Definition constraint_manager.cc:82
Class that manages the deal.II DoFHandlers.
Definition dof_manager.h:27
const std::vector< std::array< dealii::DoFHandler< dim >, 2 > > & get_dof_handlers() const
Getter function for the scalar and vector DoFHandlers.
Definition dof_manager.h:72
A little class that computes the element volume for our triangulation.
Definition invm_manager.h:28
InvMManager(const DoFManager< dim, degree > &dof_manager, const ConstraintManager< dim, degree, number > &constraint_manager, bool _calculate_scalar, bool _calculate_vector)
Constructor.
Definition invm_manager.h:39
std::vector< const SolutionVector< number > * > get_invm(const std::vector< FieldAttributes > &field_container, const std::set< unsigned int > &field_indices, unsigned int relative_level) const
Definition invm_manager.h:161
void compute_local_vector(const MatrixFree< dim, number > &_data, SolutionVector< number > &dst, const int &src, const std::pair< unsigned int, unsigned int > &cell_range) const
Definition invm_manager.h:283
static const VectorValue one
Definition invm_manager.h:348
std::vector< SolutionVector< number > > jxw_scalar
Vector that stores element volumes.
Definition invm_manager.h:341
void compute_local_scalar(const MatrixFree< dim, number > &_data, SolutionVector< number > &dst, const int &src, const std::pair< unsigned int, unsigned int > &cell_range) const
Definition invm_manager.h:264
dealii::VectorizedArray< number > ScalarValue
Definition invm_manager.h:30
unsigned int num_levels
Definition invm_manager.h:333
std::vector< SolutionVector< number > > invm_sqrt_vector
Definition invm_manager.h:346
bool calculate_vector
Definition invm_manager.h:336
std::vector< SolutionVector< number > > invm_sqrt_scalar
Definition invm_manager.h:345
std::vector< SolutionVector< number > > invm_vector
Definition invm_manager.h:344
void compute_invm()
Recompute the invm vectors.
Definition invm_manager.h:81
bool calculate_scalar
Definition invm_manager.h:335
const SolutionVector< number > & get_invm(TensorRank rank, unsigned int relative_level) const
Get the invm vector for a given rank and level.
Definition invm_manager.h:101
std::vector< SolutionVector< number > > invm_scalar
Definition invm_manager.h:343
dealii::Tensor< 1, dim, ScalarValue > VectorValue
Definition invm_manager.h:31
std::vector< SolutionVector< number > > jxw_vector
Definition invm_manager.h:342
std::vector< const SolutionVector< number > * > get_jxw(const std::vector< FieldAttributes > &field_container, const std::set< unsigned int > &field_indices, unsigned int relative_level) const
Definition invm_manager.h:176
void sqrt(SolutionVector< number > &dst, const SolutionVector< number > &src) const
Definition invm_manager.h:316
void compute_vector_invm()
Definition invm_manager.h:249
void compute_scalar_invm()
Compute element volume for the triangulation.
Definition invm_manager.h:234
const SolutionVector< number > & get_jxw(TensorRank rank, unsigned int relative_level) const
Get the integrated jxw vector for a given rank and level.
Definition invm_manager.h:123
std::vector< const SolutionVector< number > * > get_invm_sqrt(const std::vector< FieldAttributes > &field_container, const std::set< unsigned int > &field_indices, unsigned int relative_level) const
Definition invm_manager.h:191
void invert(SolutionVector< number > &dst, const SolutionVector< number > &src) const
Definition invm_manager.h:302
const SolutionVector< number > & get_invm_sqrt(TensorRank rank, unsigned int relative_level) const
Get the integrated jxw vector for a given rank and level.
Definition invm_manager.h:145
void initialize()
Initialize.
Definition invm_manager.h:211
std::vector< std::array< MatrixFree< dim, number >, 2 > > data
Matrix-free object.
Definition invm_manager.h:331
void reinit(const DoFManager< dim, degree > &dof_manager, const ConstraintManager< dim, degree, number > &constraint_manager)
Definition invm_manager.h:58
This is the main class that handles the construction and solving of user-specified PDEs.
Definition system_wide.h:24
@ Value
Use value of the variable as a criterion for refinement.
Definition grid_refiner_criterion.h:31
BlockVector< number >::BlockType SolutionVector
Typedef for solution vector.
Definition group_solution_handler.h:35
Definition conditional_ostreams.cc:20
TensorRank
Tensor rank of the field.
Definition type_enums.h:30
@ Vector
Definition type_enums.h:33
@ Scalar
Definition type_enums.h:32