7#include <deal.II/lac/la_parallel_vector.h>
8#include <deal.II/matrix_free/matrix_free.h>
10#include <prismspf/config.h>
11#include <prismspf/core/variable_attributes.h>
12#include <prismspf/user_inputs/user_input_parameters.h>
17PRISMS_PF_BEGIN_NAMESPACE
23template <
int dim,
int degree,
typename number =
double>
27 using VectorType = dealii::LinearAlgebra::distributed::Vector<number>;
28 using size_type = dealii::VectorizedArray<number>;
34 const std::map<unsigned int, variableAttributes> &_variable_attributes);
40 initialize(std::shared_ptr<dealii::MatrixFree<dim, number, size_type>> _data);
59 [[nodiscard]]
const VectorType &
60 get_invm(
const unsigned int &index)
const;
73 compute_scalar_invm();
79 compute_vector_invm();
85 const std::map<unsigned int, variableAttributes> *variable_attributes;
90 std::shared_ptr<dealii::MatrixFree<dim, number, size_type>> data;
95 VectorType invm_scalar;
100 VectorType invm_vector;
105 bool scalar_needed =
false;
110 bool vector_needed =
false;
116 unsigned int scalar_index = numbers::invalid_index;
122 unsigned int vector_index = numbers::invalid_index;
127 number tolerance = 1.0e-15;
130PRISMS_PF_END_NAMESPACE
This class handles the computation and access of the inverted mass matrix for explicit solves.
Definition invm_handler.h:25
const VectorType & get_invm(const unsigned int &index) const
Getter function for the mass matrix for the given field index (constant reference).
Definition invm_handler.cc:84
void compute_invm()
Compute the mass matrix for scalar/vector fields.
Definition invm_handler.cc:61
void recompute_invm()
Recompute the mass matrix for scalar/vector fields. This just points to compute_invm() and is used fo...
Definition invm_handler.cc:77
void initialize(std::shared_ptr< dealii::MatrixFree< dim, number, size_type > > _data)
Initialize.
Definition invm_handler.cc:48
void clear()
Clear the data so we have something that resembles the base constructor.
Definition invm_handler.cc:125