|
DUBeat 1.0.1
High-order discontinuous Galerkin methods and applications to cardiac electrophysiology
|
Class to solve the Laplace equation using the Discontinuous Galerkin method. More...
#include <laplace_dg.hpp>


Public Member Functions | |
| LaplaceDG () | |
| Constructor. | |
Public Member Functions inherited from ModelDG< basis > | |
| ModelDG (std::string model_name) | |
| Constructor. | |
| ModelDG (ModelDG< basis > &ModelDG)=default | |
| Default copy constructor. | |
| ModelDG (const ModelDG< basis > &ModelDG)=default | |
| Default const copy constructor. | |
| ModelDG (ModelDG< basis > &&ModelDG)=default | |
| Default move constructor. | |
| virtual void | declare_parameters (lifex::ParamHandler ¶ms) const override |
| Declare main parameters. | |
| virtual void | parse_parameters (lifex::ParamHandler ¶ms) override |
| Parse parameters from .prm file. | |
| unsigned int | get_dofs_per_cell () const |
| Return the number of degrees of freedom per element. | |
| virtual void | run () override |
| Run the simulation. | |
| virtual | ~ModelDG ()=default |
| Destructor. | |
Private Member Functions | |
| void | assemble_system () override |
| Assembly of the linear system. | |
Additional Inherited Members | |
Protected Member Functions inherited from ModelDG< basis > | |
| virtual void | setup_system () |
| Setup of the problem before the resolution. | |
| void | make_sparsity_pattern (const DoFHandlerDG< basis > &dof, dealii::DynamicSparsityPattern &sparsity, const dealii::AffineConstraints< double > &constraints=dealii::AffineConstraints< double >(), const bool keep_constrained_dofs=true, const dealii::types::subdomain_id subdomain_id=dealii::numbers::invalid_subdomain_id) |
| Creation of the sparsity pattern to assign to the system matrix before assembling. | |
| void | initialize_solution (lifex::LinAlg::MPI::Vector &solution_owned, lifex::LinAlg::MPI::Vector &solution) |
| To inizialize the solutions using the deal.II reinit. | |
| void | create_mesh () |
| Load the mesh from the default path. | |
| void | create_mesh (std::string mesh_path) |
| Load the mesh from a user-defined path. | |
| void | solve_system () |
| System solving. | |
| void | compute_errors (const lifex::LinAlg::MPI::Vector &solution_owned, const lifex::LinAlg::MPI::Vector &solution_ex_owned, const std::shared_ptr< dealii::Function< lifex::dim > > &u_ex, const std::shared_ptr< dealii::Function< lifex::dim > > &grad_u_ex, const char *solution_name=(char *)"u") const |
| To compute the \(L^\infty\) error, the \(L^2\) error, the \(H^1\) error and the \(DG\) error at the end of system solving, it exploits the DGComputeErrors<basis> class. | |
| virtual void | output_results (std::string output_name="solution") const |
| Output of results. | |
| virtual void | conversion_to_fem (lifex::LinAlg::MPI::Vector &sol_owned) |
| To convert a discretized solution from modal to nodal basis (does nothing if problem is already in nodal basis), in-place version. | |
| virtual lifex::LinAlg::MPI::Vector | conversion_to_fem (const lifex::LinAlg::MPI::Vector &sol_owned) const |
| To convert a discretized solution from modal to nodal basis (does nothing if problem is already in nodal basis), const version. | |
| virtual void | conversion_to_fem (lifex::LinAlg::MPI::Vector &sol_owned, const std::string fem_file_path, const unsigned int degree_fem=1, const double scaling_factor=1) const |
| virtual void | conversion_to_dub (lifex::LinAlg::MPI::Vector &sol_owned) |
| To convert a discretized solution in Dubiner basis (only for problems using Dubiner basis). | |
| virtual lifex::LinAlg::MPI::Vector | conversion_to_dub (const lifex::LinAlg::MPI::Vector &sol_owned) const |
| To convert a discretized solution in Dubiner basis (only for problems using Dubiner basis). | |
| virtual void | conversion_to_dub (lifex::LinAlg::MPI::Vector &sol_owned, const std::string fem_file_path, const unsigned int degree_fem=1, const double scaling_factor=1) const |
| void | discretize_analytical_solution (const std::shared_ptr< dealii::Function< lifex::dim > > &u_analytical, lifex::LinAlg::MPI::Vector &sol_owned) |
| Conversion of an analytical solution from FEM to basis coefficients. | |
| void | conversion_to_fem (lifex::LinAlg::MPI::Vector &sol_owned) |
| Conversion of a discretized solution from Dubiner coefficients to FEM coefficients. | |
| lifex::LinAlg::MPI::Vector | conversion_to_fem (const lifex::LinAlg::MPI::Vector &sol_owned) const |
| Conversion to FEM coefficients, const version. | |
| void | conversion_to_fem (lifex::LinAlg::MPI::Vector &sol_owned, const std::string fem_file_path, const unsigned int degree_fem, const double scaling_factor) const |
| void | conversion_to_dub (lifex::LinAlg::MPI::Vector &sol_owned) |
| Conversion of a discretized solution from FEM coefficients to Dubiner coefficients. | |
| lifex::LinAlg::MPI::Vector | conversion_to_dub (const lifex::LinAlg::MPI::Vector &sol_owned) const |
| Conversion to DUB coefficients, const version. | |
| void | conversion_to_dub (lifex::LinAlg::MPI::Vector &sol_owned, const std::string fem_file_path, const unsigned int degree_fem, const double scaling_factor) const |
| void | discretize_analytical_solution (const std::shared_ptr< dealii::Function< lifex::dim > > &u_analytical, lifex::LinAlg::MPI::Vector &sol_owned) |
| Conversion of an analytical solution from FEM to basis coefficients. | |
| void | discretize_analytical_solution (const std::shared_ptr< dealii::Function< lifex::dim > > &u_analytical, lifex::LinAlg::MPI::Vector &sol_owned) |
| Conversion of an analytical solution from FEM to basis coefficients. | |
Protected Attributes inherited from ModelDG< basis > | |
| const std::string | model_name |
| Name of the class/problem. | |
| unsigned int | prm_fe_degree |
| Polynomials degree. | |
| unsigned int | prm_n_refinements |
| Mesh refinement level (>=1). | |
| double | prm_penalty_coeff |
| DG Penalty coefficient. | |
| double | prm_stability_coeff |
| DG stabilty coefficient. | |
| double | scaling_factor = 1 |
| Scaling factor. | |
| std::shared_ptr< lifex::utils::MeshHandler > | triangulation |
| Triangulation (internal use for useful already implemented methods). | |
| unsigned int | dofs_per_cell |
| Number of degrees of freedom per cell. | |
| DoFHandlerDG< basis > | dof_handler |
| DoFHandler (internal use for useful already implemented methods). | |
| std::shared_ptr< DUBFEMHandler< basis > > | dub_fem_values |
| Member used for conversions between analytical, nodal and modal representations of the solutions. | |
| std::unique_ptr< AssembleDG< basis > > | assemble |
| Matrix assembler. | |
| lifex::utils::LinearSolverHandler< lifex::LinAlg::MPI::Vector > | linear_solver |
| Linear solver handler. | |
| lifex::utils::PreconditionerHandler | preconditioner |
| Preconditioner handler. | |
| lifex::LinAlg::MPI::SparseMatrix | matrix |
| Distributed matrix of the linear system. | |
| lifex::LinAlg::MPI::Vector | rhs |
| Distributed right hand side vector of the linear system. | |
| lifex::LinAlg::MPI::Vector | solution_owned |
| Distributed solution vector, without ghost entries. | |
| lifex::LinAlg::MPI::Vector | solution |
| Distributed solution vector, with ghost entries. | |
| lifex::LinAlg::MPI::Vector | solution_ex_owned |
| Distributed exact solution vector, without ghost entries. | |
| lifex::LinAlg::MPI::Vector | solution_ex |
| Distributed exact solution vector, without ghost entries. | |
| std::shared_ptr< lifex::utils::FunctionDirichlet > | u_ex |
| Pointer to exact solution function. | |
| std::shared_ptr< dealii::Function< lifex::dim > > | grad_u_ex |
| Pointer to exact gradient solution Function. | |
| std::shared_ptr< dealii::Function< lifex::dim > > | f_ex |
| Known forcing term. | |
| std::shared_ptr< dealii::Function< lifex::dim > > | g_n |
| Neumann boundary conditions. | |
Class to solve the Laplace equation using the Discontinuous Galerkin method.
\[ \begin{aligned} -\Delta u &= f & \quad & \text{in } \Omega = (-1, 1)^d, \\ u &= u_\mathrm{ex} & \quad & \text{on } \partial\Omega \end{aligned} \]
In particular, it can be solved using the FEM basis (basis=dealii::FE_SimplexDGP<lifex::dim>) or the Dubiner basis (basis=DUBValues<lifex::dim>).
Boundary conditions and source terms are provided assuming that the exact solution is:
\[ \begin{aligned} &d=2: \, u_\mathrm{ex}(x,y) = e^{x+y}, \hspace{6mm} (x,y) \in \Omega=(1,1)^2, \\ &d=3: \, u_\mathrm{ex}(x,y,z) = e^{x+y+z}, \hspace{6mm} (x,y,z) \in \Omega=(1,1)^3. \end{aligned} \]
Finally, \(d\) is specified in the lifex configuration.
Definition at line 161 of file laplace_dg.hpp.
|
inline |
Constructor.
Definition at line 111 of file laplace_dg.hpp.
|
overrideprivatevirtual |
Assembly of the linear system.
Implements ModelDG< basis >.
Definition at line 181 of file laplace_dg.hpp.