* mapping than a (bi-/tri-)linear one, then you need to call the functions
* <b>with</b> mapping argument should be used.
*
- * All functions take a sparse matrix object to hold the matrix to be created.
+ * All functions take a matrix object to hold the matrix to be created.
* The functions assume that the matrix is initialized with a sparsity pattern
* (SparsityPattern) corresponding to the given degree of freedom handler,
* i.e. the sparsity structure is already as needed. You can do this by
*
* See the general documentation of this namespace for more information.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
- const Mapping<dim, spacedim> &mapping,
- const DoFHandler<dim, spacedim> &dof,
- const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ const Mapping<dim, spacedim> &mapping,
+ const DoFHandler<dim, spacedim> &dof,
+ const Quadrature<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Call the create_mass_matrix() function, see above, with
* <tt>mapping=MappingQ@<dim@>(1)</tt>.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
- const DoFHandler<dim, spacedim> &dof,
- const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ const DoFHandler<dim, spacedim> &dof,
+ const Quadrature<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Assemble the @ref GlossMassMatrix "mass matrix" and a right hand side vector. If no coefficient
*
* See the general documentation of this namespace for more information.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
const Mapping<dim, spacedim> &mapping,
const DoFHandler<dim, spacedim> &dof,
const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Call the create_mass_matrix() function, see above, with
* <tt>mapping=MappingQ@<dim@>(1)</tt>.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
const DoFHandler<dim, spacedim> &dof,
const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Same function as above, but for hp-objects.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
- const hp::MappingCollection<dim, spacedim> &mapping,
- const DoFHandler<dim, spacedim> &dof,
- const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ const hp::MappingCollection<dim, spacedim> &mapping,
+ const DoFHandler<dim, spacedim> &dof,
+ const hp::QCollection<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Same function as above, but for hp-objects.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
- const DoFHandler<dim, spacedim> &dof,
- const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ const DoFHandler<dim, spacedim> &dof,
+ const hp::QCollection<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Same function as above, but for hp-objects.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
const hp::MappingCollection<dim, spacedim> &mapping,
const DoFHandler<dim, spacedim> &dof,
const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Same function as above, but for hp-objects.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
const DoFHandler<dim, spacedim> &dof,
const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
*
* See the general documentation of this namespace for more information.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
- const Mapping<dim, spacedim> &mapping,
- const DoFHandler<dim, spacedim> &dof,
- const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ const Mapping<dim, spacedim> &mapping,
+ const DoFHandler<dim, spacedim> &dof,
+ const Quadrature<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Call the create_laplace_matrix() function, see above, with
* <tt>mapping=MappingQ@<dim@>(1)</tt>.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
- const DoFHandler<dim, spacedim> &dof,
- const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ const DoFHandler<dim, spacedim> &dof,
+ const Quadrature<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Assemble the Laplace matrix and a right hand side vector. If no
*
* See the general documentation of this namespace for more information.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
const Mapping<dim, spacedim> &mapping,
const DoFHandler<dim, spacedim> &dof,
const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Call the create_laplace_matrix() function, see above, with
* <tt>mapping=MappingQ@<dim@>(1)</tt>.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
const DoFHandler<dim, spacedim> &dof,
const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Like the functions above, but for hp-objects.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
- const hp::MappingCollection<dim, spacedim> &mapping,
- const DoFHandler<dim, spacedim> &dof,
- const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ const hp::MappingCollection<dim, spacedim> &mapping,
+ const DoFHandler<dim, spacedim> &dof,
+ const hp::QCollection<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Like the functions above, but for hp-objects.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
- const DoFHandler<dim, spacedim> &dof,
- const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ const DoFHandler<dim, spacedim> &dof,
+ const hp::QCollection<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Like the functions above, but for hp-objects.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
const hp::MappingCollection<dim, spacedim> &mapping,
const DoFHandler<dim, spacedim> &dof,
const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Like the functions above, but for hp-objects.
*/
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
const DoFHandler<dim, spacedim> &dof,
const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type> *const a =
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type> *const a =
nullptr,
- const AffineConstraints<typename SparseMatrixType::value_type> &
- constraints = AffineConstraints<typename SparseMatrixType::value_type>());
+ const AffineConstraints<typename MatrixType::value_type> &constraints =
+ AffineConstraints<typename MatrixType::value_type>());
/**
* Exception
} // namespace AssemblerBoundary
- } // namespace internal
+
+
+
+ template <typename T>
+ using compress_t = decltype(std::declval<T>().compress(
+ std::declval<VectorOperation::values>()));
+
+
+
+ template <typename T>
+ constexpr bool has_compress =
+ dealii::internal::is_supported_operation<compress_t, T>;
+
+
+
+ template <
+ typename MatrixType,
+ std::enable_if_t<has_compress<MatrixType>, MatrixType> * = nullptr>
+ void
+ compress(MatrixType &A, const VectorOperation::values operation)
+ {
+ A.compress(operation);
+ }
+
+
+
+ template <
+ typename MatrixType,
+ std::enable_if_t<!has_compress<MatrixType>, MatrixType> * = nullptr>
+ void
+ compress(MatrixType &, const VectorOperation::values)
+ {
+ // nothing to do
+ }
+
+ } // namespace internal
} // namespace MatrixCreator
namespace MatrixCreator
{
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
const Mapping<dim, spacedim> &mapping,
const DoFHandler<dim, spacedim> &dof,
const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
Assert(matrix.m() == dof.n_dofs(),
ExcDimensionMismatch(matrix.m(), dof.n_dofs()));
hp::QCollection<dim> q_collection(q);
hp::MappingCollection<dim, spacedim> mapping_collection(mapping);
MatrixCreator::internal::AssemblerData::
- Scratch<dim, spacedim, typename SparseMatrixType::value_type>
+ Scratch<dim, spacedim, typename MatrixType::value_type>
assembler_data(fe_collection,
update_values | update_JxW_values |
(coefficient != nullptr ? update_quadrature_points :
mapping_collection);
MatrixCreator::internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type>
+ typename MatrixType::value_type>
copy_data;
copy_data.cell_matrix.reinit(
assembler_data.fe_collection.max_dofs_per_cell(),
dim,
spacedim,
typename DoFHandler<dim, spacedim>::active_cell_iterator,
- typename SparseMatrixType::value_type>,
- [&matrix](const internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type> &data) {
+ typename MatrixType::value_type>,
+ [&matrix](
+ const internal::AssemblerData::CopyData<typename MatrixType::value_type>
+ &data) {
MatrixCreator::internal::copy_local_to_global(
data,
&matrix,
- static_cast<Vector<typename SparseMatrixType::value_type> *>(
- nullptr));
+ static_cast<Vector<typename MatrixType::value_type> *>(nullptr));
},
assembler_data,
copy_data);
- matrix.compress(VectorOperation::add);
+ internal::compress(matrix, VectorOperation::add);
}
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
const DoFHandler<dim, spacedim> &dof,
const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
create_mass_matrix(get_default_linear_mapping(dof.get_triangulation()),
dof,
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
- const Mapping<dim, spacedim> &mapping,
- const DoFHandler<dim, spacedim> &dof,
- const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ const Mapping<dim, spacedim> &mapping,
+ const DoFHandler<dim, spacedim> &dof,
+ const Quadrature<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
Assert(matrix.m() == dof.n_dofs(),
ExcDimensionMismatch(matrix.m(), dof.n_dofs()));
hp::QCollection<dim> q_collection(q);
hp::MappingCollection<dim, spacedim> mapping_collection(mapping);
MatrixCreator::internal::AssemblerData::
- Scratch<dim, spacedim, typename SparseMatrixType::value_type>
+ Scratch<dim, spacedim, typename MatrixType::value_type>
assembler_data(fe_collection,
update_values | update_JxW_values |
update_quadrature_points,
q_collection,
mapping_collection);
MatrixCreator::internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type>
+ typename MatrixType::value_type>
copy_data;
copy_data.cell_matrix.reinit(
assembler_data.fe_collection.max_dofs_per_cell(),
dim,
spacedim,
typename DoFHandler<dim, spacedim>::active_cell_iterator,
- typename SparseMatrixType::value_type>,
- [&matrix, &rhs_vector](const internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type> &data) {
+ typename MatrixType::value_type>,
+ [&matrix, &rhs_vector](
+ const internal::AssemblerData::CopyData<typename MatrixType::value_type>
+ &data) {
MatrixCreator::internal::copy_local_to_global(data,
&matrix,
&rhs_vector);
assembler_data,
copy_data);
- matrix.compress(VectorOperation::add);
+ internal::compress(matrix, VectorOperation::add);
}
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
- const DoFHandler<dim, spacedim> &dof,
- const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ const DoFHandler<dim, spacedim> &dof,
+ const Quadrature<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
create_mass_matrix(get_default_linear_mapping(dof.get_triangulation()),
dof,
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
const hp::MappingCollection<dim, spacedim> &mapping,
const DoFHandler<dim, spacedim> &dof,
const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
Assert(matrix.m() == dof.n_dofs(),
ExcDimensionMismatch(matrix.m(), dof.n_dofs()));
ExcDimensionMismatch(matrix.n(), dof.n_dofs()));
MatrixCreator::internal::AssemblerData::
- Scratch<dim, spacedim, typename SparseMatrixType::value_type>
+ Scratch<dim, spacedim, typename MatrixType::value_type>
assembler_data(dof.get_fe_collection(),
update_values | update_JxW_values |
(coefficient != nullptr ? update_quadrature_points :
q,
mapping);
MatrixCreator::internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type>
+ typename MatrixType::value_type>
copy_data;
copy_data.cell_matrix.reinit(
assembler_data.fe_collection.max_dofs_per_cell(),
dim,
spacedim,
typename DoFHandler<dim, spacedim>::active_cell_iterator,
- typename SparseMatrixType::value_type>,
- [&matrix](const internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type> &data) {
+ typename MatrixType::value_type>,
+ [&matrix](
+ const internal::AssemblerData::CopyData<typename MatrixType::value_type>
+ &data) {
MatrixCreator::internal::copy_local_to_global(
data,
&matrix,
- static_cast<Vector<typename SparseMatrixType::value_type> *>(
- nullptr));
+ static_cast<Vector<typename MatrixType::value_type> *>(nullptr));
},
assembler_data,
copy_data);
- matrix.compress(VectorOperation::add);
+ internal::compress(matrix, VectorOperation::add);
}
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
const DoFHandler<dim, spacedim> &dof,
const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
create_mass_matrix(hp::StaticMappingQ1<dim, spacedim>::mapping_collection,
dof,
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
- const hp::MappingCollection<dim, spacedim> &mapping,
- const DoFHandler<dim, spacedim> &dof,
- const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ const hp::MappingCollection<dim, spacedim> &mapping,
+ const DoFHandler<dim, spacedim> &dof,
+ const hp::QCollection<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
Assert(matrix.m() == dof.n_dofs(),
ExcDimensionMismatch(matrix.m(), dof.n_dofs()));
ExcDimensionMismatch(matrix.n(), dof.n_dofs()));
MatrixCreator::internal::AssemblerData::
- Scratch<dim, spacedim, typename SparseMatrixType::value_type>
+ Scratch<dim, spacedim, typename MatrixType::value_type>
assembler_data(dof.get_fe_collection(),
update_values | update_JxW_values |
update_quadrature_points,
q,
mapping);
MatrixCreator::internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type>
+ typename MatrixType::value_type>
copy_data;
copy_data.cell_matrix.reinit(
assembler_data.fe_collection.max_dofs_per_cell(),
dim,
spacedim,
typename DoFHandler<dim, spacedim>::active_cell_iterator,
- typename SparseMatrixType::value_type>,
- [&matrix, &rhs_vector](const internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type> &data) {
+ typename MatrixType::value_type>,
+ [&matrix, &rhs_vector](
+ const internal::AssemblerData::CopyData<typename MatrixType::value_type>
+ &data) {
MatrixCreator::internal::copy_local_to_global(data,
&matrix,
&rhs_vector);
assembler_data,
copy_data);
- matrix.compress(VectorOperation::add);
+ internal::compress(matrix, VectorOperation::add);
}
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_mass_matrix(
- const DoFHandler<dim, spacedim> &dof,
- const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ const DoFHandler<dim, spacedim> &dof,
+ const hp::QCollection<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
create_mass_matrix(hp::StaticMappingQ1<dim, spacedim>::mapping_collection,
dof,
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
const Mapping<dim, spacedim> &mapping,
const DoFHandler<dim, spacedim> &dof,
const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
Assert(matrix.m() == dof.n_dofs(),
ExcDimensionMismatch(matrix.m(), dof.n_dofs()));
hp::QCollection<dim> q_collection(q);
hp::MappingCollection<dim, spacedim> mapping_collection(mapping);
MatrixCreator::internal::AssemblerData::
- Scratch<dim, spacedim, typename SparseMatrixType::value_type>
+ Scratch<dim, spacedim, typename MatrixType::value_type>
assembler_data(fe_collection,
update_gradients | update_JxW_values |
(coefficient != nullptr ? update_quadrature_points :
q_collection,
mapping_collection);
MatrixCreator::internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type>
+ typename MatrixType::value_type>
copy_data;
copy_data.cell_matrix.reinit(
assembler_data.fe_collection.max_dofs_per_cell(),
dim,
spacedim,
typename DoFHandler<dim, spacedim>::active_cell_iterator,
- typename SparseMatrixType::value_type>,
- [&matrix](const internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type> &data) {
+ typename MatrixType::value_type>,
+ [&matrix](
+ const internal::AssemblerData::CopyData<typename MatrixType::value_type>
+ &data) {
MatrixCreator::internal::copy_local_to_global(
data,
&matrix,
- static_cast<Vector<typename SparseMatrixType::value_type> *>(
- nullptr));
+ static_cast<Vector<typename MatrixType::value_type> *>(nullptr));
},
assembler_data,
copy_data);
- matrix.compress(VectorOperation::add);
+ internal::compress(matrix, VectorOperation::add);
}
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
const DoFHandler<dim, spacedim> &dof,
const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
create_laplace_matrix(get_default_linear_mapping(dof.get_triangulation()),
dof,
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
- const Mapping<dim, spacedim> &mapping,
- const DoFHandler<dim, spacedim> &dof,
- const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ const Mapping<dim, spacedim> &mapping,
+ const DoFHandler<dim, spacedim> &dof,
+ const Quadrature<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
Assert(matrix.m() == dof.n_dofs(),
ExcDimensionMismatch(matrix.m(), dof.n_dofs()));
hp::QCollection<dim> q_collection(q);
hp::MappingCollection<dim, spacedim> mapping_collection(mapping);
MatrixCreator::internal::AssemblerData::
- Scratch<dim, spacedim, typename SparseMatrixType::value_type>
+ Scratch<dim, spacedim, typename MatrixType::value_type>
assembler_data(fe_collection,
update_gradients | update_values | update_JxW_values |
update_quadrature_points,
q_collection,
mapping_collection);
MatrixCreator::internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type>
+ typename MatrixType::value_type>
copy_data;
copy_data.cell_matrix.reinit(
assembler_data.fe_collection.max_dofs_per_cell(),
dim,
spacedim,
typename DoFHandler<dim, spacedim>::active_cell_iterator,
- typename SparseMatrixType::value_type>,
- [&matrix, &rhs_vector](const internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type> &data) {
+ typename MatrixType::value_type>,
+ [&matrix, &rhs_vector](
+ const internal::AssemblerData::CopyData<typename MatrixType::value_type>
+ &data) {
MatrixCreator::internal::copy_local_to_global(data,
&matrix,
&rhs_vector);
assembler_data,
copy_data);
- matrix.compress(VectorOperation::add);
+ internal::compress(matrix, VectorOperation::add);
}
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
- const DoFHandler<dim, spacedim> &dof,
- const Quadrature<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ const DoFHandler<dim, spacedim> &dof,
+ const Quadrature<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
create_laplace_matrix(get_default_linear_mapping(dof.get_triangulation()),
dof,
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
const hp::MappingCollection<dim, spacedim> &mapping,
const DoFHandler<dim, spacedim> &dof,
const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
Assert(matrix.m() == dof.n_dofs(),
ExcDimensionMismatch(matrix.m(), dof.n_dofs()));
ExcDimensionMismatch(matrix.n(), dof.n_dofs()));
MatrixCreator::internal::AssemblerData::
- Scratch<dim, spacedim, typename SparseMatrixType::value_type>
+ Scratch<dim, spacedim, typename MatrixType::value_type>
assembler_data(dof.get_fe_collection(),
update_gradients | update_JxW_values |
(coefficient != nullptr ? update_quadrature_points :
q,
mapping);
MatrixCreator::internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type>
+ typename MatrixType::value_type>
copy_data;
copy_data.cell_matrix.reinit(
assembler_data.fe_collection.max_dofs_per_cell(),
dim,
spacedim,
typename DoFHandler<dim, spacedim>::active_cell_iterator,
- typename SparseMatrixType::value_type>,
- [&matrix](const internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type> &data) {
+ typename MatrixType::value_type>,
+ [&matrix](
+ const internal::AssemblerData::CopyData<typename MatrixType::value_type>
+ &data) {
MatrixCreator::internal::copy_local_to_global(
data,
&matrix,
- static_cast<Vector<typename SparseMatrixType::value_type> *>(
- nullptr));
+ static_cast<Vector<typename MatrixType::value_type> *>(nullptr));
},
assembler_data,
copy_data);
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
const DoFHandler<dim, spacedim> &dof,
const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
create_laplace_matrix(
hp::StaticMappingQ1<dim, spacedim>::mapping_collection,
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
- const hp::MappingCollection<dim, spacedim> &mapping,
- const DoFHandler<dim, spacedim> &dof,
- const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ const hp::MappingCollection<dim, spacedim> &mapping,
+ const DoFHandler<dim, spacedim> &dof,
+ const hp::QCollection<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
Assert(matrix.m() == dof.n_dofs(),
ExcDimensionMismatch(matrix.m(), dof.n_dofs()));
ExcDimensionMismatch(matrix.n(), dof.n_dofs()));
MatrixCreator::internal::AssemblerData::
- Scratch<dim, spacedim, typename SparseMatrixType::value_type>
+ Scratch<dim, spacedim, typename MatrixType::value_type>
assembler_data(dof.get_fe_collection(),
update_gradients | update_values | update_JxW_values |
update_quadrature_points,
q,
mapping);
MatrixCreator::internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type>
+ typename MatrixType::value_type>
copy_data;
copy_data.cell_matrix.reinit(
assembler_data.fe_collection.max_dofs_per_cell(),
dim,
spacedim,
typename DoFHandler<dim, spacedim>::active_cell_iterator,
- typename SparseMatrixType::value_type>,
- [&matrix, &rhs_vector](const internal::AssemblerData::CopyData<
- typename SparseMatrixType::value_type> &data) {
+ typename MatrixType::value_type>,
+ [&matrix, &rhs_vector](
+ const internal::AssemblerData::CopyData<typename MatrixType::value_type>
+ &data) {
MatrixCreator::internal::copy_local_to_global(data,
&matrix,
&rhs_vector);
- template <int dim, int spacedim, typename SparseMatrixType>
+ template <int dim, int spacedim, typename MatrixType>
void
create_laplace_matrix(
- const DoFHandler<dim, spacedim> &dof,
- const hp::QCollection<dim> &q,
- SparseMatrixType &matrix,
- const Function<spacedim, typename SparseMatrixType::value_type> &rhs,
- Vector<typename SparseMatrixType::value_type> &rhs_vector,
- const Function<spacedim, typename SparseMatrixType::value_type>
- *const coefficient,
- const AffineConstraints<typename SparseMatrixType::value_type> &constraints)
+ const DoFHandler<dim, spacedim> &dof,
+ const hp::QCollection<dim> &q,
+ MatrixType &matrix,
+ const Function<spacedim, typename MatrixType::value_type> &rhs,
+ Vector<typename MatrixType::value_type> &rhs_vector,
+ const Function<spacedim, typename MatrixType::value_type>
+ *const coefficient,
+ const AffineConstraints<typename MatrixType::value_type> &constraints)
{
create_laplace_matrix(
hp::StaticMappingQ1<dim, spacedim>::mapping_collection,