Construct a LinearProblem. More...
#include <linear_problem.h>
Public Types | |
using | ARKode = typename SUNDIALS::ARKode< typename LacType::BlockVector > |
SUNDIALS time integrator. More... | |
using | Triangulation = typename std::conditional< dim==1, parallel::shared::Triangulation< dim, spacedim >, parallel::distributed::Triangulation< dim, spacedim > >::type |
Make sure we can run also in 1d, where parallel distributed triangulations are not available, and we can only use parallel shared ones. More... | |
using | CopyData = MeshWorker::CopyData< 1, 1, 1 > |
Default CopyData object, used in the WorkStream class. More... | |
using | ScratchData = MeshWorker::ScratchData< dim, spacedim > |
Default ScratchData object, used in the workstream class. More... | |
using | BlockVectorType = typename LacType::BlockVector |
Block vector type. More... | |
using | VectorType = typename BlockVectorType::BlockType |
Vector type. More... | |
using | BlockMatrixType = typename LacType::BlockSparseMatrix |
Block matrix type. More... | |
Public Member Functions | |
LinearProblem (const std::string &component_names="u", const std::string &problem_name="") | |
Constructor. More... | |
virtual | ~LinearProblem ()=default |
Virtual destructor. More... | |
virtual void | run () |
Main entry point of the problem. More... | |
void | run_steady_state () |
Solve a steady state problem. More... | |
void | run_quasi_static () |
Solve a quasi static problem. More... | |
void | run_transient () |
Solve a dynamic problem. More... | |
virtual void | setup_transient (ARKode &arkode) |
Setup the transient problem. More... | |
virtual void | assemble_system_one_cell (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, ScratchData &scratch, CopyData ©) |
Assemble the local system matrix on cell , using scratch for FEValues and other expensive scratch objects, and store the result in the copy object. More... | |
virtual void | copy_one_cell (const CopyData ©) |
Distribute the data that has been assembled by assemble_system_on_cell() to the global matrix and rhs. More... | |
virtual void | solve () |
Solve the global system. More... | |
virtual void | estimate (Vector< float > &error_per_cell) const |
Perform a posteriori error estimation, and store the results in the error_per_cell vector. More... | |
void | mark (const Vector< float > &error_per_cell) |
According to the chosen strategy, mark some cells for refinement, possibily using the error_per_cell vector. More... | |
void | refine () |
Refine the grid. More... | |
virtual void | setup_system () |
Initial setup: distribute degrees of freedom, make all vectors and matrices of the right size, initialize functions and pointers. More... | |
virtual void | custom_estimator (Vector< float > &error_per_cell) const |
Overload this function to use a custom error estimator in the mesh refinement process. More... | |
virtual void | assemble_system () |
Actually loop over cells, and assemble the global system. More... | |
virtual void | output_results (const unsigned cycle) const |
Output the solution and the grid in a format that can be read by Paraview or Visit. More... | |
virtual void | print_system_info () const |
print some information about the current processes/mpi/ranks/etc. More... | |
![]() | |
ParameterAcceptor (const std::string §ion_name="") | |
unsigned int | get_acceptor_id () const |
virtual | ~ParameterAcceptor () override |
virtual void | declare_parameters (ParameterHandler &prm) |
virtual void | parse_parameters (ParameterHandler &prm) |
std::string | get_section_name () const |
std::vector< std::string > | get_section_path () const |
void | add_parameter (const std::string &entry, ParameterType ¶meter, const std::string &documentation="", ParameterHandler &prm_=prm, const Patterns::PatternBase &pattern=*Patterns::Tools::Convert< ParameterType >::to_pattern()) |
void | enter_subsection (const std::string &subsection) |
void | leave_subsection () |
void | enter_my_subsection (ParameterHandler &prm) |
void | leave_my_subsection (ParameterHandler &prm) |
void | serialize (Archive &ar, const unsigned int version) |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
![]() | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | serialize (Archive &ar, const unsigned int version) |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
Public Attributes | |
boost::signals2::signal< void()> | check_consistency_call_back |
Check consistency of the problem. More... | |
boost::signals2::signal< void()> | add_constraints_call_back |
A signal that is called at the end of setup_system() More... | |
boost::signals2::signal< void()> | setup_system_call_back |
A signal that is called at the end of setup_system() More... | |
boost::signals2::signal< void()> | output_results_call_back |
A signal that is called at the end of output_results() More... | |
boost::signals2::signal< void()> | assemble_system_call_back |
A signal that is called at the end of assemble_system() More... | |
boost::signals2::signal< void(ParsedTools::DataOut< dim, spacedim > &)> | add_data_vector |
Connect to this signal to add data additional vectors to the output system. More... | |
boost::signals2::signal< void(const double &time, const double &time_step, const unsigned int &time_step_number)> | advance_time_call_back |
Connect to this signal to receive time information. More... | |
const std::string | component_names |
Comma seperated names of components. More... | |
const unsigned int | n_components |
Number of components. More... | |
const std::string | problem_name |
Name of the problem to solve. More... | |
const std::string | section_name |
Name of the section to use within the parameter file. More... | |
MPI_Comm | mpi_communicator |
Global mpi communicator. More... | |
const unsigned int | mpi_rank |
The mpi rank of this process. More... | |
const unsigned int | mpi_size |
The number of mpi processes. More... | |
int | number_of_threads = 1 |
Number of threads to use for multi-threaded assembly. More... | |
unsigned int | verbosity_level = 4 |
Verbosity level of deallog. More... | |
ConditionalOStream | pcout |
Output stream, only active on process 0. More... | |
TimerOutput | timer |
Timing information. More... | |
EvolutionType | evolution_type |
Describe the type of time evolution of the problem. More... | |
ParsedTools::GridGenerator< dim, spacedim > | grid_generator |
A wrapper around GridIn, GridOut, and GridGenerator namespace. More... | |
ParsedTools::GridRefinement | grid_refinement |
Grid refinement and error estimation. More... | |
Triangulation | triangulation |
The problem triangulation. More... | |
ParsedTools::FiniteElement< dim, spacedim > | finite_element |
A wrapper around deal.II FiniteElement classes. More... | |
std::unique_ptr< Mapping< dim, spacedim > > | mapping |
The Mapping between reference and real elements. More... | |
Quadrature< dim > | cell_quadrature |
A quadrature used for cell integration. More... | |
Quadrature< dim - 1 > | face_quadrature |
A quadrature used for face integration. More... | |
DoFHandler< dim, spacedim > | dof_handler |
Handler of degrees of freedom. More... | |
AffineConstraints< double > | constraints |
Hanging nodes and essential boundary conditions. More... | |
std::vector< types::global_dof_index > | dofs_per_block |
Dofs per block. More... | |
std::vector< IndexSet > | locally_owned_dofs |
All degrees of freedom owned by this MPI process. More... | |
std::vector< IndexSet > | locally_relevant_dofs |
All degrees of freedom needed for output and error estimation. More... | |
LacType::BlockSparsityPattern | sparsity |
System sparsity pattern. More... | |
LacType::BlockSparseMatrix | matrix |
System matrix. More... | |
LacType::BlockSparseMatrix | mass_matrix |
System matrix. More... | |
LacType::BlockVector | locally_relevant_solution |
A read only copy of the solution vector used for output and error estimation. More... | |
LacType::BlockVector | solution |
Solution vector. More... | |
LacType::BlockVector | rhs |
The system right hand side. More... | |
Vector< float > | error_per_cell |
Storage for local error estimator. More... | |
ParsedLAC::InverseOperator | inverse_operator |
Inverse operator. More... | |
LacType::AMG | preconditioner |
Preconditioner. More... | |
ParsedLAC::InverseOperator | mass_inverse_operator |
Inverse operator for the mass matrix. More... | |
LacType::AMG | mass_preconditioner |
Preconditioner for the mass matrix. More... | |
ParsedTools::Function< spacedim > | forcing_term |
The actual function to use as a forcing term. More... | |
ParsedTools::Function< spacedim > | exact_solution |
The actual function to use as a exact solution when computing the errors. More... | |
ParsedTools::Function< spacedim > | initial_value |
Only used for transient problems. More... | |
ParsedTools::BoundaryConditions< spacedim > | boundary_conditions |
Boundary conditions used in this class. More... | |
ParsedTools::ConvergenceTable | error_table |
This is a wrapper around the ParsedConvergenceTable class, that allows you to specify what error to computes, and how to compute them. More... | |
ParsedTools::DataOut< dim, spacedim > | data_out |
Wrapper around the DataOut class. More... | |
double | start_time = 0.0 |
Initial time for transient and quasi stati simulations. More... | |
double | end_time = 1.0 |
Final time for transient and quasi-static simulations. More... | |
double | desired_start_step_size = .0625 |
Initial step size for transient and quasi-static simulations. More... | |
unsigned int | output_frequency = 1 |
How often to output the solution. More... | |
ParsedTools::Proxy< typename SUNDIALS::ARKode< typename LacType::BlockVector >::AdditionalData > | ark_ode_data |
Configuration used to setup transient simulations. More... | |
boost::signals2::signal< void(ARKode &)> | setup_arkode_call_back |
Signal that is triggered after creating the arkode object. More... | |
![]() | |
boost::signals2::signal< void()> | declare_parameters_call_back |
boost::signals2::signal< void()> | parse_parameters_call_back |
Static Public Attributes | |
static constexpr bool | lac_is_dealii |
True if we are using deal.II Linear Algebra Classes. More... | |
static constexpr bool | lac_is_petsc |
True if we are using PETSc Linear Algebra Classes. More... | |
static constexpr bool | lac_is_trilinos |
True if we are using Trilinos Linear Algebra Classes. More... | |
![]() | |
static ParameterHandler | prm |
Additional Inherited Members | |
![]() | |
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 void | initialize (std::istream &input_stream, ParameterHandler &prm=ParameterAcceptor::prm) |
static void | clear () |
static void | parse_all_parameters (ParameterHandler &prm=ParameterAcceptor::prm) |
static void | declare_all_parameters (ParameterHandler &prm=ParameterAcceptor::prm) |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
![]() | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
![]() | |
const std::string | section_name |
std::vector< std::string > | subsections |
Construct a LinearProblem.
Definition at line 96 of file linear_problem.h.
using PDEs::LinearProblem< dim, spacedim, LacType >::ARKode = typename SUNDIALS::ARKode<typename LacType::BlockVector> |
SUNDIALS time integrator.
Definition at line 145 of file linear_problem.h.
using PDEs::LinearProblem< dim, spacedim, LacType >::Triangulation = typename std::conditional< dim == 1, parallel::shared::Triangulation<dim, spacedim>, parallel::distributed::Triangulation<dim, spacedim> >::type |
Make sure we can run also in 1d, where parallel distributed triangulations are not available, and we can only use parallel shared ones.
Definition at line 176 of file linear_problem.h.
using PDEs::LinearProblem< dim, spacedim, LacType >::CopyData = MeshWorker::CopyData<1, 1, 1> |
Default CopyData object, used in the WorkStream class.
Definition at line 184 of file linear_problem.h.
using PDEs::LinearProblem< dim, spacedim, LacType >::ScratchData = MeshWorker::ScratchData<dim, spacedim> |
Default ScratchData object, used in the workstream class.
Definition at line 189 of file linear_problem.h.
using PDEs::LinearProblem< dim, spacedim, LacType >::BlockVectorType = typename LacType::BlockVector |
Block vector type.
Definition at line 194 of file linear_problem.h.
using PDEs::LinearProblem< dim, spacedim, LacType >::VectorType = typename BlockVectorType::BlockType |
Vector type.
Definition at line 199 of file linear_problem.h.
using PDEs::LinearProblem< dim, spacedim, LacType >::BlockMatrixType = typename LacType::BlockSparseMatrix |
Block matrix type.
Definition at line 204 of file linear_problem.h.
PDEs::LinearProblem< dim, spacedim, LacType >::LinearProblem | ( | const std::string & | component_names = "u" , |
const std::string & | problem_name = "" |
||
) |
Constructor.
Store component names and component masks.
Definition at line 38 of file linear_problem.cc.
References PDEs::LinearProblem< dim, spacedim, LacType >::custom_estimator(), MPI, and PDEs::steady_state.
|
virtualdefault |
Virtual destructor.
|
virtual |
Main entry point of the problem.
The role of this function is simply to call one of run_steady_state(), run_quasi_static() or run_transient().
Definition at line 468 of file linear_problem.cc.
References Assert, StandardExceptions::ExcNotImplemented(), PDEs::quasi_static, PDEs::steady_state, and PDEs::transient.
void PDEs::LinearProblem< dim, spacedim, LacType >::run_steady_state |
Solve a steady state problem.
Definition at line 621 of file linear_problem.cc.
References deallog, and triangulation.
void PDEs::LinearProblem< dim, spacedim, LacType >::run_quasi_static |
Solve a quasi static problem.
Definition at line 490 of file linear_problem.cc.
References DiscreteTime::advance_time(), deallog, DiscreteTime::get_current_time(), DiscreteTime::get_next_step_size(), DiscreteTime::get_previous_step_size(), DiscreteTime::get_step_number(), DiscreteTime::is_at_end(), and triangulation.
void PDEs::LinearProblem< dim, spacedim, LacType >::run_transient |
Solve a dynamic problem.
Definition at line 594 of file linear_problem.cc.
References AssertThrow, deallog, StandardExceptions::ExcMessage(), VectorTools::interpolate(), and triangulation.
|
virtual |
Setup the transient problem.
Definition at line 531 of file linear_problem.cc.
References deallog, and inverse_operator().
|
virtual |
Assemble the local system matrix on cell
, using scratch
for FEValues and other expensive scratch objects, and store the result in the copy
object.
See the documentation of WorkStream for an explanation of how to use this function.
cell | Cell on which we assemble the local matrix and rhs. |
scratch | Scratch object. |
copy | Copy object. |
Reimplemented in PDEs::MPI::Poisson< dim, spacedim >, and PDEs::LinearElasticity< dim, spacedim, LacType >.
Definition at line 249 of file linear_problem.cc.
References Assert, and StandardExceptions::ExcPureFunctionCalled().
|
virtual |
Distribute the data that has been assembled by assemble_system_on_cell() to the global matrix and rhs.
copy | The local data to distribute on the system matrix and rhs. |
Definition at line 271 of file linear_problem.cc.
|
virtual |
Solve the global system.
Reimplemented in PDEs::Stokes< dim, LacType >, PDEs::MPI::Poisson< dim, spacedim >, PDEs::LinearViscoElasticity< dim, spacedim, LacType >, and PDEs::LinearElasticity< dim, spacedim, LacType >.
Definition at line 371 of file linear_problem.cc.
References Assert, and StandardExceptions::ExcPureFunctionCalled().
|
virtual |
Perform a posteriori error estimation, and store the results in the error_per_cell
vector.
Definition at line 380 of file linear_problem.cc.
void PDEs::LinearProblem< dim, spacedim, LacType >::mark | ( | const Vector< float > & | error_per_cell | ) |
According to the chosen strategy, mark some cells for refinement, possibily using the error_per_cell
vector.
Definition at line 399 of file linear_problem.cc.
References triangulation.
void PDEs::LinearProblem< dim, spacedim, LacType >::refine |
|
virtual |
Initial setup: distribute degrees of freedom, make all vectors and matrices of the right size, initialize functions and pointers.
Definition at line 128 of file linear_problem.cc.
|
virtual |
Overload this function to use a custom error estimator in the mesh refinement process.
In order to trigger this estimator, you have to select custom
in the Error estimator
section of the parameter file.
Reimplemented in PDEs::MPI::Poisson< dim, spacedim >.
Definition at line 261 of file linear_problem.cc.
Referenced by PDEs::LinearProblem< dim, spacedim, LacType >::LinearProblem().
|
virtual |
Actually loop over cells, and assemble the global system.
Reimplemented in PDEs::LinearViscoElasticity< dim, spacedim, LacType >.
Definition at line 284 of file linear_problem.cc.
|
virtual |
Output the solution and the grid in a format that can be read by Paraview or Visit.
cycle | A number identifying the refinement cycle. |
Definition at line 421 of file linear_problem.cc.
References deallog, Utilities::int_to_string(), and Utilities::needed_digits().
|
virtual |
print some information about the current processes/mpi/ranks/etc.
Definition at line 446 of file linear_problem.cc.
References AssertThrow, deallog, StandardExceptions::ExcMessage(), MultithreadInfo::n_cores(), MultithreadInfo::n_threads(), and MultithreadInfo::set_thread_limit().
boost::signals2::signal<void()> PDEs::LinearProblem< dim, spacedim, LacType >::check_consistency_call_back |
Check consistency of the problem.
Definition at line 113 of file linear_problem.h.
|
staticconstexpr |
True if we are using deal.II Linear Algebra Classes.
Definition at line 156 of file linear_problem.h.
|
staticconstexpr |
True if we are using PETSc Linear Algebra Classes.
Definition at line 162 of file linear_problem.h.
|
staticconstexpr |
True if we are using Trilinos Linear Algebra Classes.
Definition at line 168 of file linear_problem.h.
boost::signals2::signal<void()> PDEs::LinearProblem< dim, spacedim, LacType >::add_constraints_call_back |
A signal that is called at the end of setup_system()
Definition at line 275 of file linear_problem.h.
boost::signals2::signal<void()> PDEs::LinearProblem< dim, spacedim, LacType >::setup_system_call_back |
A signal that is called at the end of setup_system()
Definition at line 280 of file linear_problem.h.
boost::signals2::signal<void()> PDEs::LinearProblem< dim, spacedim, LacType >::output_results_call_back |
A signal that is called at the end of output_results()
Definition at line 285 of file linear_problem.h.
boost::signals2::signal<void()> PDEs::LinearProblem< dim, spacedim, LacType >::assemble_system_call_back |
A signal that is called at the end of assemble_system()
Definition at line 296 of file linear_problem.h.
boost::signals2::signal<void(ParsedTools::DataOut<dim, spacedim> &)> PDEs::LinearProblem< dim, spacedim, LacType >::add_data_vector |
Connect to this signal to add data additional vectors to the output system.
Definition at line 318 of file linear_problem.h.
boost::signals2::signal<void(const double &time, const double &time_step, const unsigned int &time_step_number)> PDEs::LinearProblem< dim, spacedim, LacType >::advance_time_call_back |
Connect to this signal to receive time information.
Definition at line 326 of file linear_problem.h.
const std::string PDEs::LinearProblem< dim, spacedim, LacType >::component_names |
Comma seperated names of components.
Definition at line 332 of file linear_problem.h.
const unsigned int PDEs::LinearProblem< dim, spacedim, LacType >::n_components |
Number of components.
Definition at line 337 of file linear_problem.h.
const std::string PDEs::LinearProblem< dim, spacedim, LacType >::problem_name |
Name of the problem to solve.
Definition at line 342 of file linear_problem.h.
const std::string PDEs::LinearProblem< dim, spacedim, LacType >::section_name |
Name of the section to use within the parameter file.
Definition at line 347 of file linear_problem.h.
MPI_Comm PDEs::LinearProblem< dim, spacedim, LacType >::mpi_communicator |
Global mpi communicator.
Definition at line 352 of file linear_problem.h.
const unsigned int PDEs::LinearProblem< dim, spacedim, LacType >::mpi_rank |
The mpi rank of this process.
Definition at line 357 of file linear_problem.h.
const unsigned int PDEs::LinearProblem< dim, spacedim, LacType >::mpi_size |
The number of mpi processes.
Definition at line 362 of file linear_problem.h.
int PDEs::LinearProblem< dim, spacedim, LacType >::number_of_threads = 1 |
Number of threads to use for multi-threaded assembly.
Definition at line 367 of file linear_problem.h.
unsigned int PDEs::LinearProblem< dim, spacedim, LacType >::verbosity_level = 4 |
Verbosity level of deallog.
Definition at line 372 of file linear_problem.h.
ConditionalOStream PDEs::LinearProblem< dim, spacedim, LacType >::pcout |
Output stream, only active on process 0.
Definition at line 377 of file linear_problem.h.
|
mutable |
Timing information.
Definition at line 382 of file linear_problem.h.
EvolutionType PDEs::LinearProblem< dim, spacedim, LacType >::evolution_type |
Describe the type of time evolution of the problem.
Definition at line 387 of file linear_problem.h.
ParsedTools::GridGenerator<dim, spacedim> PDEs::LinearProblem< dim, spacedim, LacType >::grid_generator |
A wrapper around GridIn, GridOut, and GridGenerator namespace.
The action of this class is driven by the section Grid
of the parameter file :
Where you can specify what grid to generate, how to generate it, or what file to read the grid from, and to what file to write the grid to, in addition to the initial refinement of the grid.
Definition at line 409 of file linear_problem.h.
ParsedTools::GridRefinement PDEs::LinearProblem< dim, spacedim, LacType >::grid_refinement |
Grid refinement and error estimation.
This class is a wrapper around the GridRefinement namespace, and around the KellyErrorEstimator class. The action of this class is governed by the section Grid/Refinement
of the parameter file :
where you can specify the number of refinement cycles, the type of error estimator, what marking strategy to use, etc.
At the moment, local refinement is only supported on quad/hex grids. If you try to run the code with a local refinement strategy with a tria/tetra grid, an exception will be thrown at run time.
Definition at line 443 of file linear_problem.h.
Triangulation PDEs::LinearProblem< dim, spacedim, LacType >::triangulation |
The problem triangulation.
Definition at line 448 of file linear_problem.h.
ParsedTools::FiniteElement<dim, spacedim> PDEs::LinearProblem< dim, spacedim, LacType >::finite_element |
A wrapper around deal.II FiniteElement classes.
The action of this class is driven by the parameter Finite element space (u)
of the parameter file :
You should make sure that the type of finite element you specify matches the type of triangulation you are using, i.e., FE_Q is supported only on quad/hex grids, while FE_SimplexP is supported only on tri/tetra grids.
The syntax used to specify the finite element type is the same used by the FETools::get_fe_by_name() function.
Definition at line 470 of file linear_problem.h.
std::unique_ptr<Mapping<dim, spacedim> > PDEs::LinearProblem< dim, spacedim, LacType >::mapping |
The Mapping between reference and real elements.
This is a unique pointer to allow creation via parameter files.
Definition at line 477 of file linear_problem.h.
Quadrature<dim> PDEs::LinearProblem< dim, spacedim, LacType >::cell_quadrature |
A quadrature used for cell integration.
Definition at line 482 of file linear_problem.h.
Quadrature<dim - 1> PDEs::LinearProblem< dim, spacedim, LacType >::face_quadrature |
A quadrature used for face integration.
Definition at line 487 of file linear_problem.h.
DoFHandler<dim, spacedim> PDEs::LinearProblem< dim, spacedim, LacType >::dof_handler |
Handler of degrees of freedom.
Definition at line 493 of file linear_problem.h.
AffineConstraints<double> PDEs::LinearProblem< dim, spacedim, LacType >::constraints |
Hanging nodes and essential boundary conditions.
Definition at line 498 of file linear_problem.h.
std::vector<types::global_dof_index> PDEs::LinearProblem< dim, spacedim, LacType >::dofs_per_block |
Dofs per block.
Definition at line 503 of file linear_problem.h.
std::vector<IndexSet> PDEs::LinearProblem< dim, spacedim, LacType >::locally_owned_dofs |
All degrees of freedom owned by this MPI process.
Definition at line 508 of file linear_problem.h.
std::vector<IndexSet> PDEs::LinearProblem< dim, spacedim, LacType >::locally_relevant_dofs |
All degrees of freedom needed for output and error estimation.
Definition at line 513 of file linear_problem.h.
LacType::BlockSparsityPattern PDEs::LinearProblem< dim, spacedim, LacType >::sparsity |
System sparsity pattern.
Definition at line 518 of file linear_problem.h.
LacType::BlockSparseMatrix PDEs::LinearProblem< dim, spacedim, LacType >::matrix |
System matrix.
Definition at line 523 of file linear_problem.h.
LacType::BlockSparseMatrix PDEs::LinearProblem< dim, spacedim, LacType >::mass_matrix |
System matrix.
Definition at line 528 of file linear_problem.h.
LacType::BlockVector PDEs::LinearProblem< dim, spacedim, LacType >::locally_relevant_solution |
A read only copy of the solution vector used for output and error estimation.
Definition at line 534 of file linear_problem.h.
LacType::BlockVector PDEs::LinearProblem< dim, spacedim, LacType >::solution |
Solution vector.
Definition at line 539 of file linear_problem.h.
LacType::BlockVector PDEs::LinearProblem< dim, spacedim, LacType >::rhs |
The system right hand side.
Read-write vector, containing only locally owned dofs.
Definition at line 545 of file linear_problem.h.
Vector<float> PDEs::LinearProblem< dim, spacedim, LacType >::error_per_cell |
Storage for local error estimator.
This vector contains also values associated to artificial cells (i.e., it is of length triangulation.n_active_cells()
), but it is non-zero only on locally owned cells. The estimate() method only fills locally owned cells.
Definition at line 553 of file linear_problem.h.
ParsedLAC::InverseOperator PDEs::LinearProblem< dim, spacedim, LacType >::inverse_operator |
Inverse operator.
Definition at line 558 of file linear_problem.h.
LacType::AMG PDEs::LinearProblem< dim, spacedim, LacType >::preconditioner |
Preconditioner.
Definition at line 563 of file linear_problem.h.
ParsedLAC::InverseOperator PDEs::LinearProblem< dim, spacedim, LacType >::mass_inverse_operator |
Inverse operator for the mass matrix.
Definition at line 568 of file linear_problem.h.
LacType::AMG PDEs::LinearProblem< dim, spacedim, LacType >::mass_preconditioner |
Preconditioner for the mass matrix.
Definition at line 573 of file linear_problem.h.
ParsedTools::Function<spacedim> PDEs::LinearProblem< dim, spacedim, LacType >::forcing_term |
The actual function to use as a forcing term.
This is a wrapper around the ParsedFunction class, which allows you to define a function through a symbolic expression (a string) in a parameter file.
The action of this class is driven by the section Functions
, with the parameter Forcing term
:
You can use any of the numerical constants that are defined in the numbers namespace, such as PI, E, etc, as well as the constants defined at construction time in the ParsedTools::Constants class.
Definition at line 594 of file linear_problem.h.
ParsedTools::Function<spacedim> PDEs::LinearProblem< dim, spacedim, LacType >::exact_solution |
The actual function to use as a exact solution when computing the errors.
This is a wrapper around the ParsedFunction class, which allows you to define a function through a symbolic expression (a string) in a parameter file.
The action of this class is driven by the section Functions
, with the parameter Exact solution
:
You can use any of the numerical constants that are defined in the numbers namespace, such as PI, E, etc, as well as the constants defined at construction time in the ParsedTools::Constants class.
Definition at line 615 of file linear_problem.h.
ParsedTools::Function<spacedim> PDEs::LinearProblem< dim, spacedim, LacType >::initial_value |
Only used for transient problems.
Definition at line 620 of file linear_problem.h.
ParsedTools::BoundaryConditions<spacedim> PDEs::LinearProblem< dim, spacedim, LacType >::boundary_conditions |
Boundary conditions used in this class.
The action of this class is driven by the section Boundary conditions
of the parameter file :
The way ParsedTools::BoundaryConditions works in the FSI-suite is the following: for every set of boundary ids of the triangulation, you need to specify what boundary conditions are assumed to be imposed on that set. If you only want to specify one type of boundary condition (dirichlet
or neumann
) on all of the boundary, you can do so by specifying -1
as the boundary id set.
Multiple boundary conditions can be specified, but the same id should should appear only once in the parameter file (i.e., you cannot apply different types of boundary conditions on the same boundary id).
Keep in mind the following caveats:
set Boundary id sets (u) = 0, 1; 2, 3
, then boundary ids 0 and 1 will get Neumann boundary conditions, while boundary ids 2 and 3 will get Dirichlet boundary conditions.,
character, then expression for each component are separated by a semicolumn, and for each boundary id set, are separated by the %
character. For example, if you want to specify homogeneous Neumann boundary conditions, and constant Dirichlet boundary conditions you can set the following parameter: set Expressions (u) = 0 % 1
.all
, a component name, or u.n
, or u.t
to select normal component, or tangential component in a vector valued problem. For scalar problems, only the name of the component makes sense. This field allows you to control which components the given boundary condition refers to.To summarize, the following is a valid section for the example above:
This would apply Dirichlet boundary conditions on the boundary ids 2 and 3, and homogeneous Neumann boundary conditions on the boundary ids 0 and 1.
Definition at line 683 of file linear_problem.h.
|
mutable |
This is a wrapper around the ParsedConvergenceTable class, that allows you to specify what error to computes, and how to compute them.
The action of this class is driven by the section Error table
of the parameter file :
The above code, for example, would produce a convergence table that looks like
The table above can be used as-is to produce high quality pdf outputs of your error convergence rates using the file in the FSI-suite repository latex/quick_convergence_graphs/graph.tex
. For example, the above table would result in the following plot:
Definition at line 724 of file linear_problem.h.
|
mutable |
Wrapper around the DataOut class.
The action of this class is driven by the section Output
of the parameter file :
A similar structure was used in the program dof_plotter.cc.
For example, using the configuration specified above, a plot of the solution using the vtu
format would look like:
Definition at line 750 of file linear_problem.h.
double PDEs::LinearProblem< dim, spacedim, LacType >::start_time = 0.0 |
Initial time for transient and quasi stati simulations.
Definition at line 755 of file linear_problem.h.
double PDEs::LinearProblem< dim, spacedim, LacType >::end_time = 1.0 |
Final time for transient and quasi-static simulations.
Definition at line 760 of file linear_problem.h.
double PDEs::LinearProblem< dim, spacedim, LacType >::desired_start_step_size = .0625 |
Initial step size for transient and quasi-static simulations.
Definition at line 765 of file linear_problem.h.
unsigned int PDEs::LinearProblem< dim, spacedim, LacType >::output_frequency = 1 |
How often to output the solution.
Definition at line 770 of file linear_problem.h.
ParsedTools::Proxy< typename SUNDIALS::ARKode<typename LacType::BlockVector>::AdditionalData> PDEs::LinearProblem< dim, spacedim, LacType >::ark_ode_data |
Configuration used to setup transient simulations.
Definition at line 777 of file linear_problem.h.
boost::signals2::signal<void(ARKode &)> PDEs::LinearProblem< dim, spacedim, LacType >::setup_arkode_call_back |
Signal that is triggered after creating the arkode object.
Only used in transient simulations.
Definition at line 783 of file linear_problem.h.