Skip to content

Commit

Permalink
Add typename and instantiation
Browse files Browse the repository at this point in the history
  • Loading branch information
NiklasWik committed Jun 4, 2022
1 parent cc36628 commit 2d20b7d
Show file tree
Hide file tree
Showing 3 changed files with 73 additions and 47 deletions.
20 changes: 10 additions & 10 deletions include/deal.II/numerics/vector_tools_boundary.h
Original file line number Diff line number Diff line change
Expand Up @@ -718,15 +718,15 @@ namespace VectorTools
* @see
* @ref GlossBoundaryIndicator "Glossary entry on boundary indicators"
*/
template <int dim>
template <int dim, typename number, typename number2 = number>
void
project_boundary_values_div_conforming(
const DoFHandler<dim, dim> & dof_handler,
const unsigned int first_vector_component,
const Function<dim, double> &boundary_function,
const types::boundary_id boundary_component,
AffineConstraints<double> & constraints,
const Mapping<dim> & mapping);
const DoFHandler<dim, dim> & dof_handler,
const unsigned int first_vector_component,
const Function<dim, number2> &boundary_function,
const types::boundary_id boundary_component,
AffineConstraints<number> & constraints,
const Mapping<dim> & mapping);

/**
* Same as above for the hp-namespace.
Expand All @@ -736,14 +736,14 @@ namespace VectorTools
* @see
* @ref GlossBoundaryIndicator "Glossary entry on boundary indicators"
*/
template <int dim>
template <int dim, typename number, typename number2 = number>
void
project_boundary_values_div_conforming(
const DoFHandler<dim, dim> & dof_handler,
const unsigned int first_vector_component,
const Function<dim, double> & boundary_function,
const Function<dim, number2> & boundary_function,
const types::boundary_id boundary_component,
AffineConstraints<double> & constraints,
AffineConstraints<number> & constraints,
const hp::MappingCollection<dim, dim> &mapping_collection =
hp::StaticMappingQ1<dim>::mapping_collection);

Expand Down
76 changes: 43 additions & 33 deletions include/deal.II/numerics/vector_tools_boundary.templates.h
Original file line number Diff line number Diff line change
Expand Up @@ -2149,16 +2149,16 @@ namespace VectorTools
{
// This function computes the projection of the boundary function on the
// boundary in 2d.
template <typename cell_iterator>
template <typename cell_iterator, typename number, typename number2>
void
compute_face_projection_div_conforming(
const cell_iterator & cell,
const unsigned int face,
const FEFaceValues<2> & fe_values,
const unsigned int first_vector_component,
const Function<2> & boundary_function,
const Function<2, number2> & boundary_function,
const std::vector<DerivativeForm<1, 2, 2>> &jacobians,
AffineConstraints<double> & constraints)
AffineConstraints<number> & constraints)
{
// Compute the integral over the product of the normal components of
// the boundary function times the normal components of the shape
Expand All @@ -2171,9 +2171,9 @@ namespace VectorTools
1,
0,
0};
std::vector<Vector<double>> values(fe_values.n_quadrature_points,
Vector<double>(2));
Vector<double> dof_values(fe.n_dofs_per_face(face));
std::vector<Vector<number2>> values(fe_values.n_quadrature_points,
Vector<number2>(2));
Vector<number2> dof_values(fe.n_dofs_per_face(face));

// Get the values of the boundary function at the quadrature points.
{
Expand All @@ -2186,7 +2186,7 @@ namespace VectorTools
for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points;
++q_point)
{
double tmp = 0.0;
number2 tmp = 0.0;

for (unsigned int d = 0; d < 2; ++d)
tmp += normals[q_point][d] * values[q_point](d);
Expand Down Expand Up @@ -2235,32 +2235,35 @@ namespace VectorTools
}

// dummy implementation of above function for all other dimensions
template <int dim, typename cell_iterator>
template <int dim,
typename cell_iterator,
typename number,
typename number2>
void
compute_face_projection_div_conforming(
const cell_iterator &,
const unsigned int,
const FEFaceValues<dim> &,
const unsigned int,
const Function<dim> &,
const Function<dim, number2> &,
const std::vector<DerivativeForm<1, dim, dim>> &,
AffineConstraints<double> &)
AffineConstraints<number> &)
{
Assert(false, ExcNotImplemented());
}

// This function computes the projection of the boundary function on the
// boundary in 3d.
template <typename cell_iterator>
template <typename cell_iterator, typename number, typename number2>
void
compute_face_projection_div_conforming(
const cell_iterator & cell,
const unsigned int face,
const FEFaceValues<3> & fe_values,
const unsigned int first_vector_component,
const Function<3> & boundary_function,
const Function<3, number2> & boundary_function,
const std::vector<DerivativeForm<1, 3, 3>> &jacobians,
std::vector<double> & dof_values,
std::vector<number> & dof_values,
std::vector<types::global_dof_index> & projected_dofs)
{
// Compute the intergral over the product of the normal components of
Expand All @@ -2272,9 +2275,9 @@ namespace VectorTools
const unsigned int
face_coordinate_directions[GeometryInfo<3>::faces_per_cell][2] = {
{1, 2}, {1, 2}, {2, 0}, {2, 0}, {0, 1}, {0, 1}};
std::vector<Vector<double>> values(fe_values.n_quadrature_points,
Vector<double>(3));
Vector<double> dof_values_local(fe.n_dofs_per_face(face));
std::vector<Vector<number2>> values(fe_values.n_quadrature_points,
Vector<number2>(3));
Vector<number2> dof_values_local(fe.n_dofs_per_face(face));

{
const std::vector<Point<3>> &quadrature_points =
Expand All @@ -2286,7 +2289,7 @@ namespace VectorTools
for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points;
++q_point)
{
double tmp = 0.0;
number2 tmp = 0.0;

for (unsigned int d = 0; d < 3; ++d)
tmp += normals[q_point][d] * values[q_point](d);
Expand Down Expand Up @@ -2342,32 +2345,35 @@ namespace VectorTools
// dummy implementation of above
// function for all other
// dimensions
template <int dim, typename cell_iterator>
template <int dim,
typename cell_iterator,
typename number,
typename number2>
void
compute_face_projection_div_conforming(
const cell_iterator &,
const unsigned int,
const FEFaceValues<dim> &,
const unsigned int,
const Function<dim> &,
const Function<dim, number2> &,
const std::vector<DerivativeForm<1, dim, dim>> &,
std::vector<double> &,
std::vector<number> &,
std::vector<types::global_dof_index> &)
{
Assert(false, ExcNotImplemented());
}
} // namespace internals


template <int dim>
template <int dim, typename number, typename number2>
void
project_boundary_values_div_conforming(
const DoFHandler<dim> & dof_handler,
const unsigned int first_vector_component,
const Function<dim> & boundary_function,
const types::boundary_id boundary_component,
AffineConstraints<double> &constraints,
const Mapping<dim> & mapping)
const DoFHandler<dim> & dof_handler,
const unsigned int first_vector_component,
const Function<dim, number2> &boundary_function,
const types::boundary_id boundary_component,
AffineConstraints<number> & constraints,
const Mapping<dim> & mapping)
{
const unsigned int spacedim = dim;
// Interpolate the normal components
Expand Down Expand Up @@ -2430,7 +2436,9 @@ namespace VectorTools
{
AssertThrow(
dynamic_cast<const FE_RaviartThomas<dim> *>(&fe) !=
nullptr,
nullptr ||
dynamic_cast<const FE_RaviartThomasNodal<dim> *>(
&fe) != nullptr,
typename FiniteElement<
dim>::ExcInterpolationNotImplemented());
}
Expand Down Expand Up @@ -2485,7 +2493,9 @@ namespace VectorTools
{
AssertThrow(
dynamic_cast<const FE_RaviartThomas<dim> *>(&fe) !=
nullptr,
nullptr ||
dynamic_cast<const FE_RaviartThomasNodal<dim> *>(
&fe) != nullptr,
typename FiniteElement<
dim>::ExcInterpolationNotImplemented());
}
Expand Down Expand Up @@ -2529,14 +2539,14 @@ namespace VectorTools
}


template <int dim>
template <int dim, typename number, typename number2>
void
project_boundary_values_div_conforming(
const DoFHandler<dim> & dof_handler,
const unsigned int first_vector_component,
const Function<dim> & boundary_function,
const Function<dim, number2> & boundary_function,
const types::boundary_id boundary_component,
AffineConstraints<double> & constraints,
AffineConstraints<number> & constraints,
const hp::MappingCollection<dim, dim> &mapping_collection)
{
const unsigned int spacedim = dim;
Expand Down Expand Up @@ -2618,7 +2628,7 @@ namespace VectorTools
case 3:
{
const unsigned int n_dofs = dof_handler.n_dofs();
std::vector<double> dof_values(n_dofs);
std::vector<number2> dof_values(n_dofs);
std::vector<types::global_dof_index> projected_dofs(n_dofs);

for (unsigned int dof = 0; dof < n_dofs; ++dof)
Expand Down
24 changes: 20 additions & 4 deletions source/numerics/vector_tools_boundary.inst.in
Original file line number Diff line number Diff line change
Expand Up @@ -268,21 +268,37 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
AffineConstraints<std::complex<double>> &,
const hp::MappingCollection<deal_II_dimension> &);
# endif
# endif
#endif
\}
}

for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS;
number : REAL_SCALARS;
number2 : REAL_SCALARS)
{
namespace VectorTools
\{
#if deal_II_dimension == deal_II_space_dimension

# if deal_II_dimension != 1

template void
project_boundary_values_div_conforming<deal_II_dimension>(
const DoFHandler<deal_II_dimension> &,
const unsigned int,
const Function<deal_II_dimension> &,
const Function<deal_II_dimension, number2> &,
const types::boundary_id,
AffineConstraints<double> &,
AffineConstraints<number> &,
const Mapping<deal_II_dimension> &);

template void
project_boundary_values_div_conforming<deal_II_dimension>(
const DoFHandler<deal_II_dimension> &,
const unsigned int,
const Function<deal_II_dimension> &,
const Function<deal_II_dimension, number2> &,
const types::boundary_id,
AffineConstraints<double> &,
AffineConstraints<number> &,
const hp::MappingCollection<deal_II_dimension> &);
# endif
#endif
Expand Down

0 comments on commit 2d20b7d

Please sign in to comment.