Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Enable FE_DGQ::get_interpolation_matrix() for source finite element type FENothing #16496

Merged
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
4 changes: 4 additions & 0 deletions doc/news/changes/minor/20240118Schreter
Original file line number Diff line number Diff line change
@@ -0,0 +1,4 @@
Changed:
FE_DGQ::get_interpolation_matrix() now also takes as a source element FENothing.
<br>
(Magdalena Schreter, Peter Munch, 2024/01/18)
145 changes: 80 additions & 65 deletions source/fe/fe_dgq.cc
Original file line number Diff line number Diff line change
Expand Up @@ -274,78 +274,93 @@ FE_DGQ<dim, spacedim>::get_interpolation_matrix(
const FiniteElement<dim, spacedim> &x_source_fe,
FullMatrix<double> &interpolation_matrix) const
{
// this is only implemented, if the
// source FE is also a
// DGQ element
using FE = FiniteElement<dim, spacedim>;
AssertThrow((dynamic_cast<const FE_DGQ<dim, spacedim> *>(&x_source_fe) !=
nullptr),
typename FE::ExcInterpolationNotImplemented());

// ok, source is a Q element, so
// we will be able to do the work
const FE_DGQ<dim, spacedim> &source_fe =
dynamic_cast<const FE_DGQ<dim, spacedim> &>(x_source_fe);

Assert(interpolation_matrix.m() == this->n_dofs_per_cell(),
ExcDimensionMismatch(interpolation_matrix.m(),
this->n_dofs_per_cell()));
Assert(interpolation_matrix.n() == source_fe.n_dofs_per_cell(),
ExcDimensionMismatch(interpolation_matrix.n(),
source_fe.n_dofs_per_cell()));


// compute the interpolation
// matrices in much the same way as
// we do for the embedding matrices
// from mother to child.
FullMatrix<double> cell_interpolation(this->n_dofs_per_cell(),
this->n_dofs_per_cell());
FullMatrix<double> source_interpolation(this->n_dofs_per_cell(),
source_fe.n_dofs_per_cell());
FullMatrix<double> tmp(this->n_dofs_per_cell(), source_fe.n_dofs_per_cell());
for (unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
// go through the list of elements we can interpolate from
if (const FE_DGQ<dim, spacedim> *source_fe =
dynamic_cast<const FE_DGQ<dim, spacedim> *>(&x_source_fe))
{
// generate a point on this
// cell and evaluate the
// shape functions there
const Point<dim> p = this->unit_support_points[j];
for (unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
cell_interpolation(j, i) = this->poly_space->compute_value(i, p);

for (unsigned int i = 0; i < source_fe.n_dofs_per_cell(); ++i)
source_interpolation(j, i) = source_fe.poly_space->compute_value(i, p);
}
Assert(interpolation_matrix.m() == this->n_dofs_per_cell(),
ExcDimensionMismatch(interpolation_matrix.m(),
this->n_dofs_per_cell()));
Assert(interpolation_matrix.n() == source_fe->n_dofs_per_cell(),
ExcDimensionMismatch(interpolation_matrix.n(),
source_fe->n_dofs_per_cell()));


// compute the interpolation
// matrices in much the same way as
// we do for the embedding matrices
// from mother to child.
FullMatrix<double> cell_interpolation(this->n_dofs_per_cell(),
this->n_dofs_per_cell());
FullMatrix<double> source_interpolation(this->n_dofs_per_cell(),
source_fe->n_dofs_per_cell());
FullMatrix<double> tmp(this->n_dofs_per_cell(),
source_fe->n_dofs_per_cell());
for (unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
{
// generate a point on this
// cell and evaluate the
// shape functions there
const Point<dim> p = this->unit_support_points[j];
for (unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
cell_interpolation(j, i) = this->poly_space->compute_value(i, p);

for (unsigned int i = 0; i < source_fe->n_dofs_per_cell(); ++i)
source_interpolation(j, i) =
source_fe->poly_space->compute_value(i, p);
}

// then compute the
// interpolation matrix matrix
// for this coordinate
// direction
cell_interpolation.gauss_jordan();
cell_interpolation.mmult(interpolation_matrix, source_interpolation);
// then compute the
// interpolation matrix matrix
// for this coordinate
// direction
cell_interpolation.gauss_jordan();
cell_interpolation.mmult(interpolation_matrix, source_interpolation);

// cut off very small values
for (unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
for (unsigned int j = 0; j < source_fe.n_dofs_per_cell(); ++j)
if (std::fabs(interpolation_matrix(i, j)) < 1e-15)
interpolation_matrix(i, j) = 0.;
// cut off very small values
for (unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
for (unsigned int j = 0; j < source_fe->n_dofs_per_cell(); ++j)
if (std::fabs(interpolation_matrix(i, j)) < 1e-15)
interpolation_matrix(i, j) = 0.;

#ifdef DEBUG
// make sure that the row sum of
// each of the matrices is 1 at
// this point. this must be so
// since the shape functions sum up
// to 1
for (unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
{
double sum = 0.;
for (unsigned int j = 0; j < source_fe.n_dofs_per_cell(); ++j)
sum += interpolation_matrix(i, j);
// make sure that the row sum of
// each of the matrices is 1 at
// this point. this must be so
// since the shape functions sum up
// to 1
for (unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
{
double sum = 0.;
for (unsigned int j = 0; j < source_fe->n_dofs_per_cell(); ++j)
sum += interpolation_matrix(i, j);

Assert(std::fabs(sum - 1) < 5e-14 * std::max(this->degree, 1U) * dim,
ExcInternalError());
}
Assert(std::fabs(sum - 1) < 5e-14 * std::max(this->degree, 1U) * dim,
ExcInternalError());
}
#endif
}
else if (dynamic_cast<const FE_Nothing<dim> *>(&x_source_fe))
{
// the element we want to interpolate from is an FE_Nothing. this
// element represents a function that is constant zero and has no
// degrees of freedom, so the interpolation is simply a multiplication
// with a n_dofs x 0 matrix. there is nothing to do here

// we would like to verify that the number of rows and columns of
// the matrix equals this->n_dofs_per_cell() and zero. unfortunately,
// whenever we do FullMatrix::reinit(m,0), it sets both rows and
// columns to zero, instead of m and zero. thus, only test the
// number of columns
Assert(interpolation_matrix.n() == x_source_fe.n_dofs_per_cell(),
ExcDimensionMismatch(interpolation_matrix.m(),
x_source_fe.n_dofs_per_cell()));
}
else
AssertThrow(
false,
(typename FiniteElement<dim,
spacedim>::ExcInterpolationNotImplemented()));
}


Expand Down
30 changes: 30 additions & 0 deletions tests/fe/dgq_1.cc
Original file line number Diff line number Diff line change
Expand Up @@ -21,6 +21,7 @@
// result doesn't change

#include <deal.II/fe/fe_dgq.h>
#include <deal.II/fe/fe_nothing.h>
#include <deal.II/fe/fe_tools.h>

#include <string>
Expand Down Expand Up @@ -56,6 +57,31 @@ test(const unsigned int degree1, const unsigned int degree2)
}


template <int dim>
void
test_fe_nothing()
{
deallog << "FE_DGQ<" << dim << ">(1)"
<< " to FE_Nothing<" << dim << ">()" << std::endl;

FE_DGQ<dim> fe1(1);
FE_Nothing<dim> fe2;

FullMatrix<double> m;
fe1.get_interpolation_matrix(fe2, m);

for (unsigned int i = 0; i < m.m(); ++i)
{
for (unsigned int j = 0; j < m.n(); ++j)
deallog << m(i, j) << ' ';

deallog << std::endl;
}

deallog << std::endl;
}


int
main()
{
Expand All @@ -74,5 +100,9 @@ main()
for (unsigned int degree2 = 0; degree2 <= 2; ++degree2)
test<3>(degree1, degree2);

test_fe_nothing<1>();
test_fe_nothing<2>();
test_fe_nothing<3>();

return 0;
}
6 changes: 6 additions & 0 deletions tests/fe/dgq_1.output
Original file line number Diff line number Diff line change
Expand Up @@ -402,3 +402,9 @@ DEAL::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.
DEAL::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 1.00000 0.00000
DEAL::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 1.00000
DEAL::
DEAL::FE_DGQ<1>(1) to FE_Nothing<1>()
DEAL::
DEAL::FE_DGQ<2>(1) to FE_Nothing<2>()
DEAL::
DEAL::FE_DGQ<3>(1) to FE_Nothing<3>()
DEAL::
2 changes: 1 addition & 1 deletion tests/hp/solution_transfer_03.cc
Original file line number Diff line number Diff line change
Expand Up @@ -91,7 +91,7 @@ main()
data_out.write_vtu(deallog.get_file_stream());


// Interpoalte solution
// Interpolate solution
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Nice. Thank you for making these small changes too while you're there -- it is appreciated!

SolutionTransfer<2, Vector<double>> solultion_trans(dof_handler);
solultion_trans.prepare_for_coarsening_and_refinement(solution);

Expand Down