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>
19#include <prismspf/config.h>
26template <
unsigned int dim,
unsigned int degree,
typename number>
41 bool _calculate_scalar,
42 bool _calculate_vector)
43 :
num_levels(dof_manager.get_dof_handlers().size())
54 reinit(dof_manager, constraint_manager);
61 const std::vector<std::array<dealii::DoFHandler<dim>, 2>> &dof_handlers =
63 const std::vector<std::array<dealii::AffineConstraints<number>, 2>> &constraints =
65 for (
unsigned int i = 0; i <
data.size(); ++i)
67 for (
unsigned int rank = 0; rank < 2; ++rank)
70 dof_handlers[i][rank],
105 dealii::ExcInternalError(
"Requested invm that was not calculated"));
127 dealii::ExcInternalError(
"Requested invm that was not calculated"));
149 dealii::ExcInternalError(
"Requested invm that was not calculated"));
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
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)
170 out.push_back(&
get_invm(field_container[field_index].field_type, relative_level));
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
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)
185 out.push_back(&
get_jxw(field_container[field_index].field_type, relative_level));
190 std::vector<const SolutionVector<number> *>
192 const std::set<unsigned int> &field_indices,
193 unsigned int relative_level)
const
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)
201 &
get_invm_sqrt(field_container[field_index].field_type, relative_level));
236 for (
unsigned int i = 0; i <
data.size(); ++i)
251 for (
unsigned int i = 0; i <
data.size(); ++i)
266 [[maybe_unused]]
const int &src,
267 const std::pair<unsigned int, unsigned int> &cell_range)
const
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)
272 fe_eval.reinit(cell);
273 for (
unsigned int quad = 0; quad < SystemWide<dim, degree>::quadrature.size();
276 fe_eval.submit_value(dealii::make_vectorized_array<number>(1.0), quad);
278 fe_eval.integrate_scatter(dealii::EvaluationFlags::values, dst);
285 [[maybe_unused]]
const int &src,
286 const std::pair<unsigned int, unsigned int> &cell_range)
const
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)
291 fe_eval.reinit(cell);
292 for (
unsigned int quad = 0; quad < SystemWide<dim, degree>::quadrature.size();
295 fe_eval.submit_value(
one, quad);
297 fe_eval.integrate_scatter(dealii::EvaluationFlags::values, dst);
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)
309 number src_el = src.local_element(i);
311 dst.local_element(i) = src_el == 0 ? 0.0 : 1.0 / src_el;
318 src.update_ghost_values();
319 dst.update_ghost_values();
320 Assert(dst.size() == src.size(), dealii::ExcInternalError());
321 for (
unsigned int i = 0; i < src.locally_owned_size(); ++i)
323 dst.local_element(i) = std::sqrt(src.local_element(i));
330 std::vector<std::array<MatrixFree<dim, number>, 2>>
data;
350 for (
unsigned int i = 0; i < dim; ++i)
358PRISMS_PF_END_NAMESPACE
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
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:347
std::vector< SolutionVector< number > > jxw_scalar
Vector that stores element volumes.
Definition invm_manager.h:340
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:332
std::vector< SolutionVector< number > > invm_sqrt_vector
Definition invm_manager.h:345
bool calculate_vector
Definition invm_manager.h:335
std::vector< SolutionVector< number > > invm_sqrt_scalar
Definition invm_manager.h:344
std::vector< SolutionVector< number > > invm_vector
Definition invm_manager.h:343
void compute_invm()
Recompute the invm vectors.
Definition invm_manager.h:81
bool calculate_scalar
Definition invm_manager.h:334
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:342
dealii::Tensor< 1, dim, ScalarValue > VectorValue
Definition invm_manager.h:31
std::vector< SolutionVector< number > > jxw_vector
Definition invm_manager.h:341
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:330
void reinit(const DoFManager< dim, degree > &dof_manager, const ConstraintManager< dim, degree, number > &constraint_manager)
Definition invm_manager.h:58
static const dealii::QGaussLobatto< dim > quadrature
Quadrature rule.
Definition system_wide.h:41
static const dealii::MappingQ1< dim > mapping
Mappings to and from reference cell.
Definition system_wide.h:36
typename 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:52
@ Vector
Definition type_enums.h:55
@ Scalar
Definition type_enums.h:54