MoReFEM
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Protected Member Functions | Private Types | Private Member Functions | Private Attributes
MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT > Class Template Reference

CRTP base for VariationalFormulation. More...

#include <VariationalFormulation.hpp>

Inheritance diagram for MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >:
Collaboration diagram for MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >:

Public Types

using time_manager_type = TimeManagerT
 Alias to the TimeManager instance used in current model or test.
 

Public Member Functions

template<::MoReFEM::Advanced::Concept::MoReFEMDataType MoReFEMDataT>
void Init (const MoReFEMDataT &morefem_data)
 Performs the complete build of the object.
 
template<IsFactorized IsFactorizedT>
void SolveLinear (const NumberingSubset &row_numbering_subset, const NumberingSubset &col_numbering_subset, Wrappers::Petsc::print_solver_infos do_print_solver_infos=Wrappers::Petsc::print_solver_infos::yes, const std::source_location location=std::source_location::current())
 Perform a linear solve, using Petsc's KSP algorithm.
 
template<IsFactorized IsFactorizedT>
void SolveLinear (const GlobalMatrix &matrix, const GlobalVector &rhs, GlobalVector &out, Wrappers::Petsc::print_solver_infos do_print_solver_infos=Wrappers::Petsc::print_solver_infos::yes, const std::source_location location=std::source_location::current())
 Perform a linear solve, using Petsc's KSP algorithm.
 
template<enable_non_linear_solver NonLinT = NonLinearSolverT>
std::enable_if_t< NonLinT==enable_non_linear_solver::yes, void > SolveNonLinear (const NumberingSubset &row_numbering_subset, const NumberingSubset &col_numbering_subset, Wrappers::Petsc::check_convergence do_check_convergence=Wrappers::Petsc::check_convergence::yes, SNESLineSearchType line_search_type=SNESLINESEARCHBASIC, const std::source_location location=std::source_location::current())
 Perform a non-linear solve, using Petsc's Newton SNES algorithm.
 
template<VariationalFormulationNS::On AppliedOnT, BoundaryConditionMethod BoundaryConditionMethodT = BoundaryConditionMethod::pseudo_elimination>
void ApplyEssentialBoundaryCondition (const NumberingSubset &row_numbering_subset, const NumberingSubset &col_numbering_subset)
 Apply essential boundary counditions.
 
template<class LinearAlgebraT , BoundaryConditionMethod BoundaryConditionMethodT = BoundaryConditionMethod::pseudo_elimination>
void ApplyEssentialBoundaryCondition (LinearAlgebraT &linear_algebra)
 Apply essential boundary counditions on a vector or a matrix which is not a system one.
 
void WriteSolution (const TimeManagerT &time_manager, const NumberingSubset &numbering_subset) const
 Write the solution in the output directory given in the input data file.
 
const GlobalMatrixGetSystemMatrix (const NumberingSubset &row_numbering_subset, const NumberingSubset &col_numbering_subset) const
 Access to the system matrix (A in A X = R).
 
GlobalMatrixGetNonCstSystemMatrix (const NumberingSubset &row_numbering_subset, const NumberingSubset &col_numbering_subset)
 Non constant access to the system matrix (A in A X = R).
 
const GlobalVectorGetSystemRhs (const NumberingSubset &numbering_subset) const
 Access to the system rhs (R in A X = R).
 
GlobalVectorGetNonCstSystemRhs (const NumberingSubset &numbering_subset)
 Non constant access to the system rhs (R in A X = R).
 
const GlobalVectorGetSystemSolution (const NumberingSubset &numbering_subset) const
 Access to the system solution (X in A X = R).
 
GlobalVectorGetNonCstSystemSolution (const NumberingSubset &numbering_subset)
 Non constant access to the system solution (X in A X = R).
 
const TimeManagerT & GetTimeManager () const
 Object in charge of keeping track of the time-related information.
 
const FilesystemNS::DirectoryGetResultDirectory () const noexcept
 Get the directory into which all outputs are written.
 
const FilesystemNS::DirectoryGetOutputDirectory (const NumberingSubset &numbering_subset) const
 Returns the output directory into which solutions for numbering_subset are written.
 
const DirichletBoundaryCondition::vector_shared_ptrGetEssentialBoundaryConditionList () const noexcept
 List of essential boundary condition to consider.
 
const GodOfDofGetGodOfDof () const noexcept
 Get god of dof.
 
const Wrappers::MpiGetMpi () const noexcept
 

Protected Member Functions

void AllocateSystemMatrix (const NumberingSubset &row_numbering_subset, const NumberingSubset &col_numbering_subset)
 Allocate the global matrix circonscribed by the two given numbering subsets.
 
void AllocateSystemVector (const NumberingSubset &numbering_subset)
 Allocate the global vector circonscribed by the given numbering subset.
 
template<std::size_t InitialConditionIndexT, ::MoReFEM::Advanced::Concept::MoReFEMDataType MoReFEMDataT>
void SetInitialSystemSolution (const MoReFEMDataT &morefem_data, const NumberingSubset &numbering_subset, const Unknown &unknown, const FEltSpace &felt_space)
 Initialize the vector system solution with an initial condition from the input file with respect to the unknown considered.
 
template<std::size_t InitialConditionIndexT, ::MoReFEM::Advanced::Concept::MoReFEMDataType MoReFEMDataT>
void ApplyInitialConditionToVector (const MoReFEMDataT &morefem_data, const NumberingSubset &numbering_subset, const Unknown &unknown, const FEltSpace &felt_space, GlobalVector &vector)
 
const Wrappers::Petsc::SnesGetSnes () const noexcept
 Access to the Snes underlying object.
 
Wrappers::Petsc::SnesGetNonCstSnes () noexcept
 Non constant access to the Snes underlying object.
 
std::size_t GetDisplayValue () const
 Files will be written every GetDisplayValue() time iterations.
 
TimeManagerT & GetNonCstTimeManager ()
 Object in charge of keeping track of the time-related information.
 
template<enable_non_linear_solver NonLinT = NonLinearSolverT>
std::enable_if_t< NonLinT==enable_non_linear_solver::yes, const NumberingSubset & > GetRowNumberingSubset () const noexcept
 Get the NumberingSubset for the rows, used in SolveNonLinear() call.
 
template<enable_non_linear_solver NonLinT = NonLinearSolverT>
std::enable_if_t< NonLinT==enable_non_linear_solver::yes, const NumberingSubset & > GetColumnNumberingSubset () const noexcept
 Get the NumberingSubset for the columns, used in SolveNonLinear() call.
 
Special members.
template<::MoReFEM::Advanced::Concept::MoReFEMDataType MoReFEMDataT>
 VariationalFormulation (const GodOfDof &god_of_dof, DirichletBoundaryCondition::vector_shared_ptr &&boundary_condition_list, MoReFEMDataT &morefem_data)
 Constructor.
 
virtual ~VariationalFormulation ()=default
 Destructor.
 
 VariationalFormulation (const VariationalFormulation &rhs)=delete
 The copy constructor.
 
 VariationalFormulation (VariationalFormulation &&rhs)=delete
 The move constructor. - deactivated.
 
VariationalFormulationoperator= (const VariationalFormulation &rhs)=delete
 The (copy) operator=.
 
VariationalFormulationoperator= (VariationalFormulation &&rhs)=delete
 The (move) operator=.
 

Private Types

using self = VariationalFormulation<DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT>
 Alias to current class.
 

Private Member Functions

const Internal::VarfNS::GlobalMatrixStorageGetGlobalMatrixStorage () const
 Access to the class in charge of storing the global matrices for each relevant combination of numbering subsets.
 
Internal::VarfNS::GlobalMatrixStorageGetNonCstGlobalMatrixStorage ()
 Non constant access to the class in charge of storing the global matrices for each relevant combination of numbering subsets.
 
const Internal::VarfNS::GlobalVectorStorageGetRhsVectorStorage () const
 Access to the class in charge of storing the rhs vectors for each relevant numbering subset.
 
Internal::VarfNS::GlobalVectorStorageGetNonCstRhsVectorStorage ()
 Access to the class in charge of storing the rhs vectors for each relevant numbering subset.
 
const Internal::VarfNS::GlobalVectorStorageGetSolutionVectorStorage () const
 Access to the class in charge of storing the solution vectors for each relevant numbering subset.
 
Internal::VarfNS::GlobalVectorStorageGetNonCstSolutionVectorStorage ()
 Access to the class in charge of storing the solution vectors for each relevant numbering subset.
 
template<enable_non_linear_solver NonLinT = NonLinearSolverT>
std::enable_if_t< NonLinT==enable_non_linear_solver::yes, void > SetNonLinearNumberingSubset (const NumberingSubset *row, const NumberingSubset *col)
 Set both NumberingSubset for the non linear system to solve.
 
virtual void SnesMonitorFunction (PetscInt iteration_index, PetscReal norm) const
 An internal function called in SnesInterface<VariationalFormulationT>::Viewer().
 

Private Attributes

const FilesystemNS::Directoryresult_directory_
 Directory into which all outputs are written.
 
TimeManagerT & time_manager_
 Transient parameters.
 
Wrappers::Petsc::Snes::unique_ptr snes_ = nullptr
 Wrapper over Petsc solver utilities.
 
const GodOfDofgod_of_dof_
 God of dof.
 
const DirichletBoundaryCondition::vector_shared_ptr boundary_condition_list_
 List of essential boundary condition to consider.
 
std::size_t display_value_
 Files will be written every display_value_ time iterations.
 
const NumberingSubsetrow_numbering_subset_ = nullptr
 
const NumberingSubsetcolumn_numbering_subset_ = nullptr
 
const Wrappers::Mpimpi_
 Mpi object.
 
Global matrices and vectors.
Internal::VarfNS::GlobalMatrixStorage global_matrix_storage_
 List of system matrices.
 
Internal::VarfNS::GlobalVectorStorage global_rhs_storage_
 System rhs.
 
Internal::VarfNS::GlobalVectorStorage global_solution_storage_
 System solution.
 

Detailed Description

template<class DerivedT, std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
class MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >

CRTP base for VariationalFormulation.

This class is in charge of the bulk of the calculation: storage of global matrices and vectors, assembling, call to the solver.

Template Parameters
SolverIndexTIndex in the input data file of the solver to use.
DerivedTClass which will inherit from this one through CRTP. DerivedT is expected to define the following methods (otherwise the compilation will fail):
  • void SupplInit(const MoReFEMData& morefem_data) // where MoReFEMData is defined within the model instance.
  • void AllocateMatricesAndVectors()

Additionally, if NonLinearSolverT is yes the following methods are expected:

Attention
ComputeResidual() and ComputeTangent() are expected to be solely internal functions for SNESSolve(); you may use them directly if you wish but in this case, you must ensure yourself that input and output arguments have the correct layout and that output values are properly zeroed.
This class is specifically designed to be derived once (or at least that the mandatory methods are defined directly in DerivedT).

Member Typedef Documentation

◆ time_manager_type

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
using MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::time_manager_type = TimeManagerT

Alias to the TimeManager instance used in current model or test.

Constructor & Destructor Documentation

◆ VariationalFormulation() [1/3]

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
template<::MoReFEM::Advanced::Concept::MoReFEMDataType MoReFEMDataT>
MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::VariationalFormulation ( const GodOfDof & god_of_dof,
DirichletBoundaryCondition::vector_shared_ptr && boundary_condition_list,
MoReFEMDataT & morefem_data )
protected

Constructor.

The constructor by itself doesn't do much: most of the initialisation process occurs in the subsequent call to Init().

[internal] The reason of this two-step approach is a constraint from C++ language: initialisation calls methods defined in DerivedT through CRTP pattern, and this is illegal to do so in a constructor. This isn't much of an issue: VariationalFormulation is expected to be build solely within base Model class, and this class knows this call must occur (it is already coded once and for all in its internals).

Parameters
[in,out]morefem_dataThe object which encapsulates some stuff that acts as global data, such as:
  • The content of the input data.
  • Mpi related information.
  • The directory into which output is to be written.
  • Management of time iterations

The argument is noted as input and output, but most of the content is actually constant; the only part that may be modified here is the time iteration related stuff.

Parameters
[in]god_of_dofGod of dof into which the formulation works.
[in]boundary_condition_listList of Dirichlet boundary conditions to take into account in the variational formulation.

◆ ~VariationalFormulation()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
virtual MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::~VariationalFormulation ( )
protectedvirtualdefault

◆ VariationalFormulation() [2/3]

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::VariationalFormulation ( const VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT > & rhs)
protecteddelete

The copy constructor.

Parameters
[in]rhsThe object from which the construction occurs.

◆ VariationalFormulation() [3/3]

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::VariationalFormulation ( VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT > && rhs)
protecteddelete

The move constructor. - deactivated.

Parameters
[in]rhsThe object from which the construction occurs. - deactivated.

Member Function Documentation

◆ operator=() [1/2]

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
VariationalFormulation & MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::operator= ( const VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT > & rhs)
protecteddelete

The (copy) operator=.

Parameters
[in]rhsThe object from which the affectation occurs.
Returns
Reference to the object (to enable chained affectation).

◆ operator=() [2/2]

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
VariationalFormulation & MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::operator= ( VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT > && rhs)
protecteddelete

The (move) operator=.

Parameters
[in]rhsThe object from which the affectation occurs.
Returns
Reference to the object (to enable chained affectation).

◆ Init()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
template<::MoReFEM::Advanced::Concept::MoReFEMDataType MoReFEMDataT>
void MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::Init ( const MoReFEMDataT & morefem_data)

Performs the complete build of the object.

Must be called immediately after the constructor (but it is already taken care of - see constructor help for more details).

Parameters
[in]morefem_dataThe object which encapsulates some stuff that acts as global data, such as:
  • The content of the input data.
  • Mpi related information.
  • The directory into which output is to be written.
  • Management of time iterations

◆ SolveLinear() [1/2]

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
template<IsFactorized IsFactorizedT>
void MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::SolveLinear ( const NumberingSubset & row_numbering_subset,
const NumberingSubset & col_numbering_subset,
Wrappers::Petsc::print_solver_infos do_print_solver_infos = Wrappers::Petsc::print_solver_infos::yes,
const std::source_location location = std::source_location::current() )

Perform a linear solve, using Petsc's KSP algorithm.

This method does only the resolution itself; assembling or applying boundary conditions must be done beforehand.

Template Parameters
IsFactorizedTIf the same matrix has already been set up previously you should rather set this parameter to yes to avoid a call to Petsc's KSPSetOperators() and to reuse the already computed factorized matrix. However, in doubt (acceptable during development process...), the safe option is 'no'.
Parameters
[in]row_numbering_subsetSolver will be applied on the system matrix which row is given by this numbering subset.
[in]col_numbering_subsetSolver will be applied on the system matrix which column is given by this numbering subset, and with the RHS pointed out by this very same numbering subset.
Parameters
[in]locationSTL object with relevant information about the calling site (usually to help when an exception is thrown.
[in]do_print_solver_infosWhether the Petsc solver should print information about the solve in progress on standard output.

◆ SolveLinear() [2/2]

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
template<IsFactorized IsFactorizedT>
void MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::SolveLinear ( const GlobalMatrix & matrix,
const GlobalVector & rhs,
GlobalVector & out,
Wrappers::Petsc::print_solver_infos do_print_solver_infos = Wrappers::Petsc::print_solver_infos::yes,
const std::source_location location = std::source_location::current() )

Perform a linear solve, using Petsc's KSP algorithm.

This method does only the resolution itself; assembling or applying boundary conditions must be done beforehand.

Template Parameters
IsFactorizedTIf the same matrix has already been set up previously you should rather set this parameter to yes to avoid a call to Petsc's KSPSetOperators() and to reuse the already computed factorized matrix. However, in doubt (acceptable during development process...), the safe option is 'no'.

This version is broader that the overload above, as you can give your own matrix and vector rather than using the so-called 'system' ones.

Parameters
[in]matrixMatrix part of the equation to solve (e.g. A in 'A X = R').
[in]rhsRhs part of the equation to solve (e.g. R in 'A X = R').
[out]outSolution part of the equation to solve (e.g. X in 'A X = R'). It must have been properly allocated; only the content is filled.
Parameters
[in]locationSTL object with relevant information about the calling site (usually to help when an exception is thrown.
[in]do_print_solver_infosWhether the Petsc solver should print information about the solve in progress on standard output.

◆ SolveNonLinear()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
template<enable_non_linear_solver NonLinT = NonLinearSolverT>
std::enable_if_t< NonLinT==enable_non_linear_solver::yes, void > MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::SolveNonLinear ( const NumberingSubset & row_numbering_subset,
const NumberingSubset & col_numbering_subset,
Wrappers::Petsc::check_convergence do_check_convergence = Wrappers::Petsc::check_convergence::yes,
SNESLineSearchType line_search_type = SNESLINESEARCHBASIC,
const std::source_location location = std::source_location::current() )

Perform a non-linear solve, using Petsc's Newton SNES algorithm.

Parameters
[in]row_numbering_subsetSolver will be applied on the system matrix which row is given by this numbering subset.
[in]col_numbering_subsetSolver will be applied on the system matrix which column is given by this numbering subset, and with the RHS pointed out by this very same numbering subset.
Parameters
[in]locationSTL object with relevant information about the calling site (usually to help when an exception is thrown.
[in]do_check_convergenceCheck the convergence if needed. If it does not converge throw an error.
[in]line_search_typeThe line search strategy to use. Up to PETSc 3.12, I didn't interact at all with the underlying line search method, but this version introduced a change that made the hyperelastic model converge with difficulty for the static case - the reason was that for some reason the backtrace strategy was used in place of the basic one. Default is set to the basic method; the available methods are on PETSc online documentation: https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESLineSearchSetType.html

◆ ApplyEssentialBoundaryCondition() [1/2]

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
template<VariationalFormulationNS::On AppliedOnT, BoundaryConditionMethod BoundaryConditionMethodT = BoundaryConditionMethod::pseudo_elimination>
void MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::ApplyEssentialBoundaryCondition ( const NumberingSubset & row_numbering_subset,
const NumberingSubset & col_numbering_subset )

Apply essential boundary counditions.

[internal] Natural boundary conditions are managed by dedicated operators (their 'naturalness' is fully used!)

Template Parameters
AppliedOnTWhether the essential boundary conditions are applied on matrix, vector or both.
BoundaryConditionMethodChoice of the method used to apply Dirichlet boundary conditions.
Parameters
[in]row_numbering_subsetSolver will be applied on the system matrix which row is given by this numbering subset.
[in]col_numbering_subsetSolver will be applied on the system matrix which column is given by this numbering subset, and with the RHS pointed out by this very same numbering subset.

◆ ApplyEssentialBoundaryCondition() [2/2]

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
template<class LinearAlgebraT , BoundaryConditionMethod BoundaryConditionMethodT = BoundaryConditionMethod::pseudo_elimination>
void MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::ApplyEssentialBoundaryCondition ( LinearAlgebraT & linear_algebra)

Apply essential boundary counditions on a vector or a matrix which is not a system one.

[internal] Natural boundary conditions are managed by dedicated operators (their 'naturalness' is fully used!)

Template Parameters
BoundaryConditionMethodChoice of the method used to apply Dirichlet boundary conditions.
LinearAlgebraTEither GlobalMatrix or GlobalVector.
Parameters
[in]linear_algebraVector or matrix to apply the boundary conditions on.

◆ WriteSolution()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
void MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::WriteSolution ( const TimeManagerT & time_manager,
const NumberingSubset & numbering_subset ) const

Write the solution in the output directory given in the input data file.

Parameters
[in,out]time_managerObject in charge of keeping track of the time-related information.
[in]numbering_subsetNumbering subset onto which solution is described.

◆ GetSystemMatrix()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
const GlobalMatrix & MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::GetSystemMatrix ( const NumberingSubset & row_numbering_subset,
const NumberingSubset & col_numbering_subset ) const

Access to the system matrix (A in A X = R).

Parameters
[in]row_numbering_subsetMatrix required is the one which row is described by this numbering subset.
[in]col_numbering_subsetMatrix required is the one which column is described by this numbering subset.
Returns
Reference to the adequate system matrix.

◆ GetNonCstSystemMatrix()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
GlobalMatrix & MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::GetNonCstSystemMatrix ( const NumberingSubset & row_numbering_subset,
const NumberingSubset & col_numbering_subset )

Non constant access to the system matrix (A in A X = R).

Parameters
[in]row_numbering_subsetMatrix required is the one which row is described by this numbering subset.
[in]col_numbering_subsetMatrix required is the one which column is described by this numbering subset.
Returns
Reference to the adequate system matrix.

◆ GetSystemRhs()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
const GlobalVector & MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::GetSystemRhs ( const NumberingSubset & numbering_subset) const

Access to the system rhs (R in A X = R).

Parameters
[in]numbering_subsetVector required is the one described by this numbering subset.
Returns
Reference to the required vector.

◆ GetNonCstSystemRhs()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
GlobalVector & MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::GetNonCstSystemRhs ( const NumberingSubset & numbering_subset)

Non constant access to the system rhs (R in A X = R).

Parameters
[in]numbering_subsetVector required is the one described by this numbering subset.
Returns
Reference to the required vector.

◆ GetSystemSolution()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
const GlobalVector & MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::GetSystemSolution ( const NumberingSubset & numbering_subset) const

Access to the system solution (X in A X = R).

Parameters
[in]numbering_subsetVector required is the one described by this numbering subset.
Returns
Reference to the required vector.

◆ GetNonCstSystemSolution()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
GlobalVector & MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::GetNonCstSystemSolution ( const NumberingSubset & numbering_subset)

Non constant access to the system solution (X in A X = R).

Parameters
[in]numbering_subsetVector required is the one described by this numbering subset.
Returns
Reference to the required vector.

◆ GetTimeManager()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
const TimeManagerT & MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::GetTimeManager ( ) const

Object in charge of keeping track of the time-related information.

Returns
Reference to the TimeManager.

◆ GetResultDirectory()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
const FilesystemNS::Directory & MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::GetResultDirectory ( ) const
noexcept

Get the directory into which all outputs are written.

It is assumed to already exist when this call is made.

Returns
Directory into which all outputs are written.

◆ GetOutputDirectory()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
const FilesystemNS::Directory & MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::GetOutputDirectory ( const NumberingSubset & numbering_subset) const

Returns the output directory into which solutions for numbering_subset are written.

Parameters
[in]numbering_subsetNumberingSubset of the solution for which the output directory is requested.
Returns
Path to the output directory into which solutions related to numbering_subset are written.

◆ AllocateSystemMatrix()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
void MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::AllocateSystemMatrix ( const NumberingSubset & row_numbering_subset,
const NumberingSubset & col_numbering_subset )
protected

Allocate the global matrix circonscribed by the two given numbering subsets.

This method is merely a wrapper over MoReFEM::AllocateGlobalMatrix() specialized for the system matrices within the VariationalFormulation.

Parameters
[in]row_numbering_subsetMatrix required is the one which row is described by this numbering subset.
[in]col_numbering_subsetMatrix required is the one which column is described by this numbering subset.

This is done once and for all during the initialisation phase of MoReFEM; no others should be allocated in the core of the calculation.

The data have already been reduced to processor-wise when this operation is performed (or for that matter when VariationalFormulation are expected to be built, i.e. typically in the method SupplInitialize() to define in each model instance).

The pattern of the matrix is expected to have been computed in the GodOfDof before the reduction of data.

◆ AllocateSystemVector()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
void MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::AllocateSystemVector ( const NumberingSubset & numbering_subset)
protected

Allocate the global vector circonscribed by the given numbering subset.

This method is merely a wrapper over MoReFEM::AllocateGlobalVector() specialized for the system vectors within the VariationalFormulation.

Parameters
[in]numbering_subsetVector required is the one described by this numbering subset.

This is done once and for all during the initialisation phase of MoReFEM; no others should be allocated in the core of the calculation.

The data have already been reduced to processor-wise when this operation is performed (or for that matter when VariationalFormulation are expected to be built, i.e. typically in method SupplInitialize() of the model).

[internal] Currently both rhs and solutions are allocated by this call.

◆ SetInitialSystemSolution()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
template<std::size_t InitialConditionIndexT, ::MoReFEM::Advanced::Concept::MoReFEMDataType MoReFEMDataT>
void MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::SetInitialSystemSolution ( const MoReFEMDataT & morefem_data,
const NumberingSubset & numbering_subset,
const Unknown & unknown,
const FEltSpace & felt_space )
protected

Initialize the vector system solution with an initial condition from the input file with respect to the unknown considered.

Boundary conditions related to numbering_subset are also applied by this method.

Todo
#1728 Couldn't this functionality be fulfilled with a source operator?
Template Parameters
InitialConditionIndexTIndex used to specify which inittial condition to consider from the input data file.
InputDataTType of the structure that holds the model settings.
Parameters
[in]morefem_dataThe object which encapsulates some stuff that acts as global data, such as:
  • The content of the input data.
  • Mpi related information.
  • The directory into which output is to be written.
  • Management of time iterations
[in]numbering_subsetNumberingSubset related to the unknown in the felt_space.
[in]felt_spaceFEltSpace onto which the vector is to be initialized.
[in]unknownUnknown considered in the vector.

◆ ApplyInitialConditionToVector()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
template<std::size_t InitialConditionIndexT, ::MoReFEM::Advanced::Concept::MoReFEMDataType MoReFEMDataT>
void MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::ApplyInitialConditionToVector ( const MoReFEMDataT & morefem_data,
const NumberingSubset & numbering_subset,
const Unknown & unknown,
const FEltSpace & felt_space,
GlobalVector & vector )
protected

Todo
#1728 Couldn't this functionality be fulfilled with a source operator?
Template Parameters
InitialConditionIndexTIndex used to specify which inittial condition to consider from the input data file.
InputDataTType of the structure that holds the model settings.
Parameters
[in]morefem_dataThe object which encapsulates some stuff that acts as global data, such as:
  • The content of the input data.
  • Mpi related information.
  • The directory into which output is to be written.
  • Management of time iterations
[in]numbering_subsetNumberingSubset related to the unknown in the felt_space.
[in]felt_spaceFEltSpace onto which the vector is to be initialized.
[in]unknownUnknown considered in the vector.

Fill the content of vector with data from an InitialCondition object.

Parameters
[in,out]vectorVector which value must be initialized. This vector must already been properly allocated, and itss associated numbering subset must match numbering_subset.

◆ GetNonCstTimeManager()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
TimeManagerT & MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::GetNonCstTimeManager ( )
protected

Object in charge of keeping track of the time-related information.

Returns
Reference to the TimeManager.

◆ SetNonLinearNumberingSubset()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
template<enable_non_linear_solver NonLinT = NonLinearSolverT>
std::enable_if_t< NonLinT==enable_non_linear_solver::yes, void > MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::SetNonLinearNumberingSubset ( const NumberingSubset * row,
const NumberingSubset * col )
private

Set both NumberingSubset for the non linear system to solve.

Parameters
[in]rowThe NumberingSubset used to index the rows of the system matrix.
[in]colThe NumberingSubset used to index the columns of the system matrix.

◆ SnesMonitorFunction()

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
virtual void MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::SnesMonitorFunction ( PetscInt iteration_index,
PetscReal norm ) const
privatevirtual

An internal function called in SnesInterface<VariationalFormulationT>::Viewer().

This Viewer function is what is provided to SNESMonitorSet.

Current function describes the required text to be printed at each Newton iteration. It is virtual as a model develop may want to substitute its own text to the one provided here by default; if so he must in its VariationalFormulation class which inherits from this one declare:

virtual void SnesMonitorFunction(PetscInt iteration_index, PetscReal norm) const override;
virtual void SnesMonitorFunction(PetscInt iteration_index, PetscReal norm) const
An internal function called in SnesInterface<VariationalFormulationT>::Viewer().
Parameters
[in]iteration_indexIndex of the Newton iteration; this is the quantity provided directly by PETSc (hence the type).
[in]normNorm of the current residual; this is the quantity provided directly by PETSc (hence the type).

As this function is virtual I can't restrict it to the case on linear solver is enabled at compile time; however an assert and/or an exit is emitted at runtime if it is nonetheless used.

◆ GetMpi()

template<class DerivedT >
const Wrappers::Mpi & MoReFEM::Crtp::CrtpMpi< DerivedT >::GetMpi ( ) const
noexceptinherited

Read-only access to underlying Mpi object.

Returns
Constant accessor to underlying Mpi object.

Field Documentation

◆ snes_

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
Wrappers::Petsc::Snes::unique_ptr MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::snes_ = nullptr
private

Wrapper over Petsc solver utilities.

[internal] For conveniency same class is used for all cases (including mere KSP solve).

◆ row_numbering_subset_

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
const NumberingSubset* MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::row_numbering_subset_ = nullptr
private

The NumberingSubset for the rows, used in SolveNonLinear() call. Not used if non linear solver is not enabled.

◆ column_numbering_subset_

template<class DerivedT , std::size_t SolverIndexT, TIME_MANAGER_TEMPLATE_KEYWORD TimeManagerT, enable_non_linear_solver NonLinearSolverT = enable_non_linear_solver::no>
const NumberingSubset* MoReFEM::VariationalFormulation< DerivedT, SolverIndexT, TimeManagerT, NonLinearSolverT >::column_numbering_subset_ = nullptr
private

The NumberingSubset for the columns, used in SolveNonLinear() call. Not used if non linear solver is not enabled.


The documentation for this class was generated from the following file: