Theory and Practice of FEM -- a parametrized Poisson Solver
poisson.cc
Go to the documentation of this file.
1 /* ---------------------------------------------------------------------
2  *
3  * Copyright (C) 1999 - 2020 by the deal.II authors
4  *
5  * This file is part of the deal.II library.
6  *
7  * The deal.II library is free software; you can use it, redistribute
8  * it, and/or modify it under the terms of the GNU Lesser General
9  * Public License as published by the Free Software Foundation; either
10  * version 2.1 of the License, or (at your option) any later version.
11  * The full text of the license can be found in the file LICENSE.md at
12  * the top level directory of deal.II.
13  *
14  * ---------------------------------------------------------------------
15  *
16  * Authors: Wolfgang Bangerth, 1999,
17  * Guido Kanschat, 2011
18  * Luca Heltai, 2021
19  */
20 #include "poisson.h"
21 
25 
27 
30 
33 
35 
36 using namespace dealii;
37 
38 template <int dim>
40  : timer(std::cout, TimerOutput::summary, TimerOutput::cpu_and_wall_times)
41  , dof_handler(triangulation)
42  , solver_control("Solver control", 1000, 1e-12, 1e-12)
43 {
44  TimerOutput::Scope timer_section(timer, "constructor");
45  add_parameter("Finite element degree", fe_degree);
46  add_parameter("Mapping degree", mapping_degree);
47  add_parameter("Number of global refinements", n_refinements);
48  add_parameter("Output filename", output_filename);
49  add_parameter("Forcing term expression", forcing_term_expression);
50  add_parameter("Dirichlet boundary condition expression",
52  add_parameter("Coefficient expression", coefficient_expression);
53  add_parameter("Number of threads", number_of_threads);
54  add_parameter("Exact solution expression", exact_solution_expression);
55  add_parameter("Neumann boundary condition expression",
57 
58  add_parameter("Local pre-refinement grid size expression",
60 
61  add_parameter("Dirichlet boundary ids", dirichlet_ids);
62  add_parameter("Neumann boundary ids", neumann_ids);
63 
64  add_parameter("Problem constants", constants);
65  add_parameter("Grid generator function", grid_generator_function);
66  add_parameter("Grid generator arguments", grid_generator_arguments);
67  add_parameter("Number of refinement cycles", n_refinement_cycles);
68 
69  add_parameter("Estimator type",
71  "",
72  this->prm,
73  Patterns::Selection("exact|kelly|residual"));
74 
75  add_parameter("Marking strategy",
77  "",
78  this->prm,
79  Patterns::Selection("global|fixed_fraction|fixed_number"));
80 
81  add_parameter("Coarsening and refinement factors",
83 
84  add_parameter("Use direct solver", use_direct_solver);
85 
86  this->prm.enter_subsection("Error table");
89 }
90 
91 
92 template <int dim>
93 void
94 Poisson<dim>::initialize(const std::string &filename)
95 {
96  TimerOutput::Scope timer_section(timer, "initialize");
98  "last_used_parameters.prm",
100 }
101 
102 
103 
104 template <int dim>
105 void
106 Poisson<dim>::parse_string(const std::string &parameters)
107 {
108  TimerOutput::Scope timer_section(timer, "parse_string");
111 }
112 
113 
114 
115 template <int dim>
116 void
118 {
119  TimerOutput::Scope timer_section(timer, "make_grid");
120 
121  const auto vars = dim == 1 ? "x" : dim == 2 ? "x,y" : "x,y,z";
122  pre_refinement.initialize(vars, pre_refinement_expression, constants);
124  grid_generator_function,
125  grid_generator_arguments);
126 
127  for (unsigned int i = 0; i < n_refinements; ++i)
128  {
129  for (const auto &cell : triangulation.active_cell_iterators())
130  if (pre_refinement.value(cell->center()) < cell->diameter())
131  cell->set_refine_flag();
132  triangulation.execute_coarsening_and_refinement();
133  }
134 
135  std::cout << "Number of active cells: " << triangulation.n_active_cells()
136  << std::endl;
137 }
138 
139 
140 
141 template <int dim>
142 void
144 {
145  TimerOutput::Scope timer_section(timer, "refine_grid");
146  // Cells have been marked in the mark() method.
147  triangulation.execute_coarsening_and_refinement();
148 }
149 
150 
151 
152 template <int dim>
153 void
155 {
156  TimerOutput::Scope timer_section(timer, "setup_system");
157  if (!fe)
158  {
159  fe = std::make_unique<FE_Q<dim>>(fe_degree);
160  mapping = std::make_unique<MappingQGeneric<dim>>(mapping_degree);
161  const auto vars = dim == 1 ? "x" : dim == 2 ? "x,y" : "x,y,z";
162  forcing_term.initialize(vars, forcing_term_expression, constants);
163  coefficient.initialize(vars, coefficient_expression, constants);
164  exact_solution.initialize(vars, exact_solution_expression, constants);
165 
166  dirichlet_boundary_condition.initialize(
167  vars, dirichlet_boundary_conditions_expression, constants);
168 
169  neumann_boundary_condition.initialize(
170  vars, neumann_boundary_conditions_expression, constants);
171  }
172 
173  dof_handler.distribute_dofs(*fe);
174  std::cout << "Number of degrees of freedom: " << dof_handler.n_dofs()
175  << std::endl;
176  constraints.clear();
177  DoFTools::make_hanging_node_constraints(dof_handler, constraints);
178 
179  for (const auto &id : dirichlet_ids)
181  *mapping, dof_handler, id, dirichlet_boundary_condition, constraints);
182  constraints.close();
183 
184 
185  DynamicSparsityPattern dsp(dof_handler.n_dofs());
186  DoFTools::make_sparsity_pattern(dof_handler, dsp, constraints, false);
187  sparsity_pattern.copy_from(dsp);
188  system_matrix.reinit(sparsity_pattern);
189  solution.reinit(dof_handler.n_dofs());
190  system_rhs.reinit(dof_handler.n_dofs());
191  error_per_cell.reinit(triangulation.n_active_cells());
192 }
193 
194 
195 
196 template <int dim>
197 void
199  const typename DoFHandler<dim>::active_cell_iterator &cell,
200  ScratchData & scratch,
201  CopyData & copy)
202 {
203  auto &cell_matrix = copy.matrices[0];
204  auto &cell_rhs = copy.vectors[0];
205 
206  cell->get_dof_indices(copy.local_dof_indices[0]);
207 
208  const auto &fe_values = scratch.reinit(cell);
209  cell_matrix = 0;
210  cell_rhs = 0;
211 
212  for (const unsigned int q_index : fe_values.quadrature_point_indices())
213  {
214  for (const unsigned int i : fe_values.dof_indices())
215  for (const unsigned int j : fe_values.dof_indices())
216  cell_matrix(i, j) +=
217  (coefficient.value(fe_values.quadrature_point(q_index)) * // a(x_q)
218  fe_values.shape_grad(i, q_index) * // grad phi_i(x_q)
219  fe_values.shape_grad(j, q_index) * // grad phi_j(x_q)
220  fe_values.JxW(q_index)); // dx
221  for (const unsigned int i : fe_values.dof_indices())
222  cell_rhs(i) +=
223  (fe_values.shape_value(i, q_index) * // phi_i(x_q)
224  forcing_term.value(fe_values.quadrature_point(q_index)) * // f(x_q)
225  fe_values.JxW(q_index)); // dx
226  }
227 
228  if (cell->at_boundary())
229  // for(const auto face: cell->face_indices())
230  for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
231  if (neumann_ids.find(cell->face(f)->boundary_id()) != neumann_ids.end())
232  {
233  auto &fe_face_values = scratch.reinit(cell, f);
234  for (const unsigned int q_index :
235  fe_face_values.quadrature_point_indices())
236  for (const unsigned int i : fe_face_values.dof_indices())
237  cell_rhs(i) += fe_face_values.shape_value(i, q_index) *
238  neumann_boundary_condition.value(
239  fe_face_values.quadrature_point(q_index)) *
240  fe_face_values.JxW(q_index);
241  }
242 }
243 
244 
245 
246 template <int dim>
247 void
249 {
250  constraints.distribute_local_to_global(copy.matrices[0],
251  copy.vectors[0],
252  copy.local_dof_indices[0],
253  system_matrix,
254  system_rhs);
255 }
256 
257 
258 
259 template <int dim>
260 void
262  const typename DoFHandler<dim>::active_cell_iterator &begin,
263  const typename DoFHandler<dim>::active_cell_iterator &end)
264 {
265  QGauss<dim> quadrature_formula(fe->degree + 1);
266  QGauss<dim - 1> face_quadrature_formula(fe->degree + 1);
267 
268  ScratchData scratch(*mapping,
269  *fe,
270  quadrature_formula,
273  face_quadrature_formula,
276 
277  CopyData copy(fe->n_dofs_per_cell());
278 
279  static Threads::Mutex assemble_mutex;
280 
281  for (auto cell = begin; cell != end; ++cell)
282  {
283  assemble_system_one_cell(cell, scratch, copy);
284  assemble_mutex.lock();
285  copy_one_cell(copy);
286  assemble_mutex.unlock();
287  }
288 }
289 
290 
291 
292 template <int dim>
293 void
295 {
296  TimerOutput::Scope timer_section(timer, "assemble_system");
297  QGauss<dim> quadrature_formula(fe->degree + 1);
298  QGauss<dim - 1> face_quadrature_formula(fe->degree + 1);
299 
300  ScratchData scratch(*mapping,
301  *fe,
302  quadrature_formula,
305  face_quadrature_formula,
308 
309  CopyData copy(fe->n_dofs_per_cell());
310 
311  WorkStream::run(dof_handler.begin_active(),
312  dof_handler.end(),
313  *this,
316  scratch,
317  copy);
318 }
319 
320 
321 
322 template <int dim>
323 void
325 {
326  TimerOutput::Scope timer_section(timer, "solve");
327  if (use_direct_solver == true)
328  {
329  SparseDirectUMFPACK system_matrix_inverse;
330  system_matrix_inverse.initialize(system_matrix);
331  system_matrix_inverse.vmult(solution, system_rhs);
332  }
333  else
334  {
335  SolverCG<Vector<double>> solver(solver_control);
336 #ifdef DEAL_II_WITH_TRILINOS
338  amg.initialize(system_matrix);
339  solver.solve(system_matrix, solution, system_rhs, amg);
340 #else
341  solver.solve(system_matrix, solution, system_rhs, PreconditionIdentity());
342 #endif
343  }
344  constraints.distribute(solution);
345 }
346 
347 
348 
349 template <int dim>
350 void
352 {
353  TimerOutput::Scope timer_section(timer, "estimate");
354  if (estimator_type == "exact")
355  {
356  error_per_cell = 0;
357  QGauss<dim> quad(fe->degree + 1);
359  dof_handler,
360  solution,
361  exact_solution,
362  error_per_cell,
363  quad,
365  }
366  else if (estimator_type == "kelly")
367  {
368  std::map<types::boundary_id, const Function<dim> *> neumann;
369  for (const auto id : neumann_ids)
370  neumann[id] = &neumann_boundary_condition;
371 
372  QGauss<dim - 1> face_quad(fe->degree + 1);
374  dof_handler,
375  face_quad,
376  neumann,
377  solution,
378  error_per_cell,
379  ComponentMask(),
380  &coefficient);
381  }
382  else if (estimator_type == "residual")
383  {
384  // h_T || f+\Delta u_h ||_0,T
385  // + \sum over faces
386  // 1/2 (h_F)^{1/2} || [n.\nabla u_h] ||_0,F
387 
388  QGauss<dim - 1> face_quad(fe->degree + 1);
389  QGauss<dim> quad(fe->degree + 1);
390 
391  std::map<types::boundary_id, const Function<dim> *> neumann;
392  for (const auto id : neumann_ids)
393  neumann[id] = &neumann_boundary_condition;
394 
396  dof_handler,
397  face_quad,
398  neumann,
399  solution,
400  error_per_cell,
401  ComponentMask(),
402  &coefficient);
403 
404  FEValues<dim> fe_values(*mapping,
405  *fe,
406  quad,
409 
410  std::vector<double> local_laplacians(quad.size());
411 
412 
413  double residual_L2_norm = 0;
414 
415  unsigned int cell_index = 0;
416  for (const auto &cell : dof_handler.active_cell_iterators())
417  {
418  fe_values.reinit(cell);
419 
420  fe_values.get_function_laplacians(solution, local_laplacians);
421  residual_L2_norm = 0;
422  for (const auto q_index : fe_values.quadrature_point_indices())
423  {
424  const auto arg =
425  (local_laplacians[q_index] +
426  forcing_term.value(fe_values.quadrature_point(q_index)));
427  residual_L2_norm += arg * arg * fe_values.JxW(q_index);
428  }
429  error_per_cell[cell_index] +=
430  cell->diameter() * std::sqrt(residual_L2_norm);
431 
432  ++cell_index;
433  }
434  }
435  else
436  {
437  AssertThrow(false, ExcNotImplemented());
438  }
439  auto global_estimator = error_per_cell.l2_norm();
440  error_table.add_extra_column("estimator", [global_estimator]() {
441  return global_estimator;
442  });
443  error_table.error_from_exact(*mapping, dof_handler, solution, exact_solution);
444 }
445 
446 
447 
448 template <int dim>
449 void
451 {
452  TimerOutput::Scope timer_section(timer, "mark");
453  if (marking_strategy == "global")
454  {
455  for (const auto &cell : triangulation.active_cell_iterators())
456  cell->set_refine_flag();
457  }
458  else if (marking_strategy == "fixed_fraction")
459  {
462  error_per_cell,
463  coarsening_and_refinement_factors.second,
464  coarsening_and_refinement_factors.first);
465  }
466  else if (marking_strategy == "fixed_number")
467  {
470  error_per_cell,
471  coarsening_and_refinement_factors.second,
472  coarsening_and_refinement_factors.first);
473  }
474  else
475  {
476  Assert(false, ExcInternalError());
477  }
478 }
479 
480 template <int dim>
481 void
482 Poisson<dim>::output_results(const unsigned cycle) const
483 {
484  TimerOutput::Scope timer_section(timer, "output_results");
485  DataOut<dim> data_out;
486  DataOutBase::VtkFlags flags;
487  flags.write_higher_order_cells = true;
488  data_out.set_flags(flags);
489  data_out.attach_dof_handler(dof_handler);
490  data_out.add_data_vector(solution, "solution");
491  auto interpolated_exact = solution;
492  VectorTools::interpolate(*mapping,
493  dof_handler,
494  exact_solution,
495  interpolated_exact);
496  data_out.add_data_vector(interpolated_exact, "exact");
497  data_out.add_data_vector(error_per_cell, "estimator");
498  data_out.build_patches(*mapping,
499  std::max(mapping_degree, fe_degree),
501  std::string fname = output_filename + "_" + std::to_string(cycle) + ".vtu";
502  std::ofstream output(fname);
503  data_out.write_vtu(output);
504 }
505 
506 
507 
508 template <int dim>
509 void
511 {
512  if (number_of_threads != -1 && number_of_threads > 0)
514  static_cast<unsigned int>(number_of_threads));
515 
516  std::cout << "Number of cores : " << MultithreadInfo::n_cores() << std::endl
517  << "Number of threads: " << MultithreadInfo::n_threads()
518  << std::endl;
519 }
520 
521 
522 
523 template <int dim>
524 void
526 {
527  print_system_info();
528  make_grid();
529  for (unsigned int cycle = 0; cycle < n_refinement_cycles; ++cycle)
530  {
531  setup_system();
532  assemble_system();
533  solve();
534  estimate();
535  output_results(cycle);
536  if (cycle < n_refinement_cycles - 1)
537  {
538  mark();
539  refine_grid();
540  }
541  }
542  error_table.output_table(std::cout);
543 }
544 
545 template class Poisson<1>;
546 template class Poisson<2>;
547 template class Poisson<3>;
void write_vtu(std::ostream &out) const
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
virtual void build_patches(const unsigned int n_subdivisions=0)
void get_function_laplacians(const InputVector &fe_function, std::vector< typename InputVector::value_type > &laplacians) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const
const Point< spacedim > & quadrature_point(const unsigned int q) const
double JxW(const unsigned int quadrature_point) const
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, dim >, level_dof_access >> &cell)
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandlerType &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, typename InputVector::value_type > * > &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
const FEValues< dim, spacedim > & reinit(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
static unsigned int n_cores()
static unsigned int n_threads()
static void set_thread_limit(const unsigned int max_threads=numbers::invalid_unsigned_int)
static void initialize(const std::string &filename="", const std::string &output_filename="", const ParameterHandler::OutputStyle output_style_for_output_filename=ParameterHandler::Short, ParameterHandler &prm=ParameterAcceptor::prm, const ParameterHandler::OutputStyle output_style_for_filename=ParameterHandler::DefaultStyle)
static ParameterHandler prm
static void parse_all_parameters(ParameterHandler &prm=ParameterAcceptor::prm)
void add_parameter(const std::string &entry, ParameterType &parameter, const std::string &documentation="", ParameterHandler &prm_=prm, const Patterns::PatternBase &pattern=*Patterns::Tools::Convert< ParameterType >::to_pattern())
void leave_subsection()
virtual void parse_input_from_string(const std::string &s, const std::string &last_line="", const bool skip_undefined=false)
void enter_subsection(const std::string &subsection)
void add_parameters(ParameterHandler &prm)
unsigned int size() const
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
void initialize(const SparsityPattern &sparsity_pattern)
void vmult(Vector< double > &dst, const Vector< double > &src) const
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
Poisson()
Definition: poisson.cc:39
std::set< types::boundary_id > neumann_ids
Definition: poisson.h:167
void estimate()
Definition: poisson.cc:351
std::set< types::boundary_id > dirichlet_ids
Definition: poisson.h:166
unsigned int fe_degree
Definition: poisson.h:159
void assemble_system_one_cell(const typename DoFHandler< dim >::active_cell_iterator &cell, ScratchData &scratch, CopyData &copy)
Definition: poisson.cc:198
std::map< std::string, double > constants
Definition: poisson.h:175
std::string neumann_boundary_conditions_expression
Definition: poisson.h:173
void assemble_system_on_range(const typename DoFHandler< dim >::active_cell_iterator &begin, const typename DoFHandler< dim >::active_cell_iterator &end)
Definition: poisson.cc:261
std::string exact_solution_expression
Definition: poisson.h:171
void refine_grid()
Definition: poisson.cc:143
void solve()
Definition: poisson.cc:324
void setup_system()
Definition: poisson.cc:154
std::string forcing_term_expression
Definition: poisson.h:169
std::string pre_refinement_expression
Definition: poisson.h:174
void parse_string(const std::string &par)
Definition: poisson.cc:106
void output_results(const unsigned cycle) const
Definition: poisson.cc:482
void mark()
Definition: poisson.cc:450
void initialize(const std::string &filename)
Definition: poisson.cc:94
int number_of_threads
Definition: poisson.h:164
std::pair< double, double > coarsening_and_refinement_factors
Definition: poisson.h:148
ParsedConvergenceTable error_table
Definition: poisson.h:180
std::string marking_strategy
Definition: poisson.h:147
std::string grid_generator_function
Definition: poisson.h:177
std::string coefficient_expression
Definition: poisson.h:170
bool use_direct_solver
Definition: poisson.h:182
unsigned int mapping_degree
Definition: poisson.h:160
unsigned int n_refinement_cycles
Definition: poisson.h:162
void print_system_info()
Definition: poisson.cc:510
std::string grid_generator_arguments
Definition: poisson.h:178
void run()
Definition: poisson.cc:525
std::string output_filename
Definition: poisson.h:163
std::string dirichlet_boundary_conditions_expression
Definition: poisson.h:172
void copy_one_cell(const CopyData &copy)
Definition: poisson.cc:248
TimerOutput timer
Definition: poisson.h:133
void assemble_system()
Definition: poisson.cc:294
void make_grid()
Definition: poisson.cc:117
std::string estimator_type
Definition: poisson.h:146
unsigned int n_refinements
Definition: poisson.h:161
static ::ExceptionBase & ExcInternalError()
#define Assert(cond, exc)
static ::ExceptionBase & ExcNotImplemented()
#define AssertThrow(cond, exc)
void make_sparsity_pattern(const DoFHandlerType &dof_handler, SparsityPatternType &sparsity_pattern, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
void make_hanging_node_constraints(const DoFHandlerType &dof_handler, AffineConstraints< number > &constraints)
update_hessians
update_values
update_JxW_values
update_gradients
update_quadrature_points
void generate_from_name_and_arguments(Triangulation< dim, spacedim > &tria, const std::string &grid_generator_function_name, const std::string &grid_generator_function_arguments)
void refine_and_coarsen_fixed_number(Triangulation< dim, spacedim > &triangulation, const Vector< Number > &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
void refine_and_coarsen_fixed_fraction(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double top_fraction, const double bottom_fraction, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double >> &velocity, const double factor=1.)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
VectorType::value_type * begin(VectorType &V)
VectorType::value_type * end(VectorType &V)
void interpolate(const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const ComponentMask &component_mask=ComponentMask())
void interpolate_boundary_values(const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask=ComponentMask())
void integrate_difference(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim, typename InVector::value_type > &exact_solution, OutVector &difference, const Quadrature< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=nullptr, const double exponent=2.)
void run(const std::vector< std::vector< Iterator >> &colored_iterators, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
void copy(const T *begin, const T *end, U *dest)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation