* ---------------------------------------------------------------------
*/
-// @sect3{Include files}
-
-#include <deal.II/base/convergence_table.h>
-#include <deal.II/base/tensor.h>
+// @sect3{This test considers two elements which share a common, regular face.
+// The logic of make_flux_sparsity_pattern() loops over all faces and filters
+// some of them according to a predicate face_has_flux_coupling(). We check
+// that the only internal face in the current setup is visited exactly once and
+// that the predicate face_has_flux_coupling() is evaluated once for this face.}
#include <deal.II/dofs/dof_tools.h>
#include <deal.II/fe/fe_q.h>
-#include <deal.II/fe/fe_values.h>
-#include <deal.II/fe/mapping_q1.h>
-#include <deal.II/grid/filtered_iterator.h>
#include <deal.II/grid/grid_generator.h>
-#include <deal.II/grid/grid_out.h>
-#include <deal.II/grid/tria.h>
-#include <deal.II/lac/affine_constraints.h>
#include <deal.II/lac/dynamic_sparsity_pattern.h>
-#include <deal.II/lac/solver_control.h>
-#include <deal.II/lac/sparse_direct.h>
-#include <deal.II/lac/sparsity_pattern.h>
-
-#include <deal.II/numerics/data_out.h>
-#include <deal.II/numerics/vector_tools.h>
-
-#include <fstream>
#include "../tests.h"
+using namespace dealii;
+template <int dim>
+Triangulation<dim>
+make_two_elements()
+{
+ Triangulation<dim> triangulation;
+ GridGenerator::hyper_cube(triangulation, -2.0, 2.0);
+ triangulation.begin_active()->set_refine_flag(
+ RefinementCase<dim>::cut_axis(0));
+ triangulation.execute_coarsening_and_refinement();
+ return triangulation;
+}
-// @sect3{The main() function}
-// test that a face is visited once
-using namespace dealii;
template <int dim>
bool
-is_face_on_OY(const typename Triangulation<dim>::active_cell_iterator &cell,
- const unsigned int face_index)
+is_face_on_OY(const typename DoFHandler<dim>::active_cell_iterator &cell,
+ const unsigned int face_index)
{
- bool on_OY = (std::abs(cell->face(face_index)->center()(0)) < 0.01);
- deallog.get_file_stream()
- << "This sentence should appear once when the corresponding face is visited only once"
- << std::endl;
- deallog.get_file_stream() << cell->index() << std::endl;
- return on_OY;
+ deallog
+ << "This sentence should appear once when the corresponding face is visited only once on cell "
+ << cell->index() << std::endl;
+ return (std::abs(cell->face(face_index)->center()(0)) < 0.01);
}
template <int dim>
void
-make_two_elements(Triangulation<dim> &triangulation)
+create_and_output_flux_sparsity_with_filter(
+ DoFHandler<dim> & dof_handler,
+ std::function<bool(const typename DoFHandler<dim>::active_cell_iterator,
+ const unsigned int)> filter)
{
- GridGenerator::hyper_cube(triangulation, -2.0, 2.0);
- triangulation.begin_active()->set_refine_flag(
- RefinementCase<dim>::cut_axis(0));
- triangulation.execute_coarsening_and_refinement();
+ DynamicSparsityPattern dsp(dof_handler.n_dofs(), dof_handler.n_dofs());
+ DoFTools::make_flux_sparsity_pattern(
+ dof_handler,
+ dsp,
+ AffineConstraints<double>(),
+ true,
+ Table<2, DoFTools::Coupling>(1, 1, new auto(DoFTools::always)),
+ Table<2, DoFTools::Coupling>(1, 1, new auto(DoFTools::always)),
+ numbers::invalid_subdomain_id,
+ [&](const auto &cell, const unsigned int face_index) {
+ return filter(cell, face_index);
+ });
+ dsp.print(deallog.get_file_stream());
}
template <int dim>
void
-create_flux_pattern_with_filter(DoFHandler<dim> & dof_handler,
- DynamicSparsityPattern &dsp)
+check()
{
- Table<2, DoFTools::Coupling> cell_coupling(1, 1);
- Table<2, DoFTools::Coupling> face_coupling(1, 1);
- cell_coupling[0][0] = DoFTools::always;
- face_coupling[0][0] = DoFTools::always;
-
- const AffineConstraints<double> constraints;
- const bool keep_constrained_dofs = true;
-
- const auto face_has_flux_coupling = [&](const auto & cell,
- const unsigned int face_index) {
- return is_face_on_OY<dim>(cell, face_index);
- };
-
- DoFTools::make_flux_sparsity_pattern(dof_handler,
- dsp,
- constraints,
- keep_constrained_dofs,
- cell_coupling,
- face_coupling,
- numbers::invalid_subdomain_id,
- face_has_flux_coupling);
+ // create a square with two elements
+ const Triangulation<dim> tria = make_two_elements<dim>();
+ // Generate Q1 dofs for the mesh grid
+ DoFHandler<dim> dof_handler(tria);
+ dof_handler.distribute_dofs(FE_Q<dim>(1));
+ // create and output the sparsity pattern by using a filter
+ create_and_output_flux_sparsity_with_filter<dim>(dof_handler,
+ is_face_on_OY<dim>);
}
-
int
main()
{
initlog();
- const int dim = 2;
-
- // create a square with two elements
- Triangulation<dim> triangulation;
- make_two_elements(triangulation);
-
- deallog << "dealii::Sparsity" << std::endl;
- {
- // Generate Q1 dofs for the mesh grid
- DoFHandler<dim> dof_handler(triangulation);
- const FE_Q<dim> finite_element(1);
- dof_handler.distribute_dofs(finite_element);
-
- // Compute the sparsity pattern specifying a face filter
- DynamicSparsityPattern dsp(dof_handler.n_dofs(), dof_handler.n_dofs());
- create_flux_pattern_with_filter(dof_handler, dsp);
- }
+ deallog.push("2d");
+ check<2>();
+ deallog.pop();
+ deallog.push("3d");
+ check<3>();
+ deallog.pop();
}
-DEAL::dealii::Sparsity
-This sentence should appear once when the corresponding face is visited only once
-1
+DEAL:2d::This sentence should appear once when the corresponding face is visited only once on cell 1
+[0,0,1,2,3,4,5]
+[1,0,1,2,3,4,5]
+[2,0,1,2,3,4,5]
+[3,0,1,2,3,4,5]
+[4,0,1,2,3,4,5]
+[5,0,1,2,3,4,5]
+DEAL:3d::This sentence should appear once when the corresponding face is visited only once on cell 1
+[0,0,1,2,3,4,5,6,7,8,9,10,11]
+[1,0,1,2,3,4,5,6,7,8,9,10,11]
+[2,0,1,2,3,4,5,6,7,8,9,10,11]
+[3,0,1,2,3,4,5,6,7,8,9,10,11]
+[4,0,1,2,3,4,5,6,7,8,9,10,11]
+[5,0,1,2,3,4,5,6,7,8,9,10,11]
+[6,0,1,2,3,4,5,6,7,8,9,10,11]
+[7,0,1,2,3,4,5,6,7,8,9,10,11]
+[8,0,1,2,3,4,5,6,7,8,9,10,11]
+[9,0,1,2,3,4,5,6,7,8,9,10,11]
+[10,0,1,2,3,4,5,6,7,8,9,10,11]
+[11,0,1,2,3,4,5,6,7,8,9,10,11]
\ No newline at end of file