deal.II version GIT relicensing-2167-g9622207b8f 2024-11-21 12:40:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
parsed_convergence_table.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2019 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_base_parsed_convergence_table_h
16#define dealii_base_parsed_convergence_table_h
17
18#include <deal.II/base/config.h>
19
24
26
27#include <deal.II/fe/mapping.h>
28
29#include <deal.II/grid/tria.h>
30
31#include <deal.II/lac/vector.h>
32
34
36
37#ifndef DOXYGEN
39#endif
40
124{
125public:
163 const std::vector<std::string> &component_names = {"u"},
164 const std::vector<std::set<VectorTools::NormType>> &list_of_error_norms = {
166
232 const std::vector<std::string> &component_names,
233 const std::vector<std::set<VectorTools::NormType>> &list_of_error_norms,
234 const double exponent,
235 const std::set<std::string> &extra_columns,
236 const std::string &rate_key,
237 const std::string &rate_mode,
238 const std::string &error_file_name,
239 const unsigned int precision,
240 const bool compute_error);
241
247 void
249
261 template <int dim, int spacedim, typename VectorType>
262 void
264 const VectorType &solution,
265 const Function<spacedim> &exact,
266 const Function<spacedim> *weight = nullptr);
267
271 template <int dim, int spacedim, typename VectorType>
272 void
274 const DoFHandler<dim, spacedim> &vspace,
275 const VectorType &solution,
276 const Function<spacedim> &exact,
277 const Function<spacedim> *weight = nullptr);
278
345 void
346 add_extra_column(const std::string &column_name,
347 const std::function<double()> &custom_function,
348 const bool compute_rate = true);
349
353 template <int dim, int spacedim, typename VectorType>
354 void
356 const VectorType &,
357 const VectorType &,
358 const Function<spacedim> *weight = nullptr);
359
363 template <int dim, int spacedim, typename VectorType>
364 void
367 const VectorType &,
368 const VectorType &,
369 const Function<spacedim> *weight = nullptr);
370
376 void
377 output_table(std::ostream &out);
378
385 void
386 output_table();
387
388private:
392 void
394
398 const std::vector<std::string> component_names;
399
403 const std::vector<std::string> unique_component_names;
404
408 const std::vector<ComponentMask> unique_component_masks;
409
413 std::map<std::string, std::pair<std::function<double()>, bool>>
415
419 std::vector<std::set<VectorTools::NormType>> norms_per_unique_component;
420
424 double exponent;
425
430
434 std::set<std::string> extra_columns;
435
439 std::string rate_key;
440
444 std::string rate_mode;
445
449 unsigned int precision;
450
454 std::string error_file_name;
455
461};
462
463
464
465#ifndef DOXYGEN
466// ============================================================
467// Template functions
468// ============================================================
469template <int dim, int spacedim, typename VectorType>
470void
472 const VectorType &solution1,
473 const VectorType &solution2,
474 const Function<spacedim> *weight)
475{
476 AssertThrow(solution1.size() == solution2.size(),
477 ExcDimensionMismatch(solution1.size(), solution2.size()));
478 VectorType solution(solution1);
479 solution -= solution2;
481 dh,
482 solution,
484 weight);
485}
486
487
488
489template <int dim, int spacedim, typename VectorType>
490void
493 const VectorType &solution1,
494 const VectorType &solution2,
495 const Function<spacedim> *weight)
496{
497 AssertThrow(solution1.size() == solution2.size(),
498 ExcDimensionMismatch(solution1.size(), solution2.size()));
499 VectorType solution(solution1);
500 solution -= solution2;
502 mapping,
503 dh,
504 solution,
506 weight);
507}
508
509
510
511template <int dim, int spacedim, typename VectorType>
512void
514 const VectorType &solution,
515 const Function<spacedim> &exact,
516 const Function<spacedim> *weight)
517{
519 dh,
520 solution,
521 exact,
522 weight);
523}
524
525
526
527template <int dim, int spacedim, typename VectorType>
528void
531 const VectorType &solution,
532 const Function<spacedim> &exact,
533 const Function<spacedim> *weight)
534{
535 const auto n_components = component_names.size();
536
537 if (compute_error)
538 {
539 AssertDimension(exact.n_components, n_components);
540 AssertDimension(dh.get_fe().n_components(), n_components);
541
544 const unsigned int n_dofs = dh.n_dofs();
545
546 for (const auto &col : extra_columns)
547 if (col == "cells")
548 {
549 table.add_value("cells", n_active_cells);
550 table.set_tex_caption("cells", "\\# cells");
551 table.set_tex_format("cells", "r");
552 }
553 else if (col == "dofs")
554 {
555 table.add_value("dofs", n_dofs);
556 table.set_tex_caption("dofs", "\\# dofs");
557 table.set_tex_format("dofs", "r");
558 }
559
560 // A vector of zero std::functions with n_components components
561 const std::vector<std::function<double(const Point<spacedim> &)>>
562 zero_components(n_components,
563 [](const Point<spacedim> &) { return 0.0; });
564
565 // The default weight function, with n_components components
566 std::vector<std::function<double(const Point<spacedim> &)>>
567 weight_components(n_components,
568 [](const Point<spacedim> &) { return 1.0; });
569
570 if (weight != nullptr)
571 {
572 if (weight->n_components == 1)
573 {
574 for (auto &f : weight_components)
575 f = [&](const Point<spacedim> &p) { return weight->value(p); };
576 }
577 else
578 {
579 AssertDimension(weight->n_components, n_components);
580 for (unsigned int i = 0; i < n_components; ++i)
581 weight_components[i] = [&](const Point<spacedim> &p) {
582 return weight->value(p, i);
583 };
584 }
585 }
586
587 for (unsigned int i = 0; i < norms_per_unique_component.size(); ++i)
588 {
589 std::map<VectorTools::NormType, double> errors;
590
591 const auto &norms = norms_per_unique_component[i];
592 const auto &mask = unique_component_masks[i];
593
594 // Simple case first
595 if (norms.empty())
596 continue;
597
598 auto components_expr = zero_components;
599 for (unsigned int j = 0; j < n_components; ++j)
600 if (mask[j] == true)
601 components_expr[j] = weight_components[j];
602
604 components_expr);
605
606 Vector<float> difference_per_cell(
608
609 QGauss<dim> q_gauss((dh.get_fe().degree + 1) * 2);
610
611 for (const auto &norm : norms)
612 {
613 difference_per_cell = 0;
615 dh,
616 solution,
617 exact,
618 difference_per_cell,
619 q_gauss,
620 norm,
621 &select_component,
622 exponent);
623
625 dh.get_triangulation(), difference_per_cell, norm, exponent);
626
627 std::string name = unique_component_names[i] + "_" +
629 std::string latex_name = "@f$\\| " + unique_component_names[i] +
630 " - " + unique_component_names[i] +
631 "_h \\|_{" +
632 Patterns::Tools::to_string(norm) + "}@f$";
633
634 table.add_value(name, errors[norm]);
636 table.set_scientific(name, true);
637 table.set_tex_caption(name, latex_name);
638 }
639 }
640
641 for (const auto &extra_col : extra_column_functions)
642 {
643 const double custom_error = extra_col.second.first();
644
645 std::string name = extra_col.first;
646 table.add_value(name, custom_error);
648 table.set_scientific(name, true);
649 }
650 }
651}
652
653#endif
654
656
657#endif
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const Triangulation< dim, spacedim > & get_triangulation() const
types::global_dof_index n_dofs() const
const unsigned int degree
Definition fe_data.h:452
unsigned int n_components() const
const unsigned int n_components
Definition function.h:163
virtual RangeNumberType value(const Point< dim > &p, const unsigned int component=0) const
Abstract base class for mapping classes.
Definition mapping.h:318
The ParsedConvergenceTable class.
void add_extra_column(const std::string &column_name, const std::function< double()> &custom_function, const bool compute_rate=true)
void difference(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &, const VectorType &, const VectorType &, const Function< spacedim > *weight=nullptr)
void difference(const DoFHandler< dim, spacedim > &, const VectorType &, const VectorType &, const Function< spacedim > *weight=nullptr)
void add_parameters(ParameterHandler &prm)
void error_from_exact(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &vspace, const VectorType &solution, const Function< spacedim > &exact, const Function< spacedim > *weight=nullptr)
void error_from_exact(const DoFHandler< dim, spacedim > &vspace, const VectorType &solution, const Function< spacedim > &exact, const Function< spacedim > *weight=nullptr)
const std::vector< std::string > component_names
std::set< std::string > extra_columns
const std::vector< std::string > unique_component_names
std::vector< std::set< VectorTools::NormType > > norms_per_unique_component
const std::vector< ComponentMask > unique_component_masks
std::map< std::string, std::pair< std::function< double()>, bool > > extra_column_functions
Definition point.h:111
void set_tex_format(const std::string &key, const std::string &format="c")
void add_value(const std::string &key, const T value)
void set_tex_caption(const std::string &key, const std::string &tex_caption)
void set_scientific(const std::string &key, const bool scientific)
void set_precision(const std::string &key, const unsigned int precision)
virtual types::global_cell_index n_global_active_cells() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define AssertThrow(cond, exc)
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition mapping.cc:294
Tpetra::Vector< Number, LO, GO, NodeType< MemorySpace > > VectorType
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
Definition divergence.h:471
std::string to_string(const T &t)
Definition patterns.h:2451
double compute_global_error(const Triangulation< dim, spacedim > &tria, const InVector &cellwise_error, const NormType &norm, const double exponent=2.)
void integrate_difference(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const ReadVector< Number > &fe_function, const Function< spacedim, Number > &exact_solution, OutVector &difference, const Quadrature< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=nullptr, const double exponent=2.)
unsigned int n_active_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
Definition tria.cc:14878