|
| v0.14.0
|
Go to the documentation of this file.
17 : cOre(const_cast<
Core &>(core)) {}
53 auto copy_dm_struture = [&](
auto simple_dm) {
64 dm = copy_dm_struture(
simple->getDM());
69 auto set_dm_section = [&](
auto dm) {
73 CHKERR DMSetSection(dm, section);
78 boost::shared_ptr<FEMethod>
null;
113 auto copy_dm_struture = [&](
auto simple_dm) {
124 dm = copy_dm_struture(
simple->getDM());
129 auto set_dm_section = [&](
auto dm) {
133 CHKERR DMSetSection(dm, section);
136 CHKERR set_dm_section(dm);
140 boost::shared_ptr<FEMethod>
null;
165 CHKERR SNESSetDM(snes, dm);
186 "TS solver handling not implemented");
196 auto copy_dm_struture = [&](
auto simple_dm) {
207 dm = copy_dm_struture(
simple->getDM());
212 auto set_dm_section = [&](
auto dm) {
216 CHKERR DMSetSection(dm, section);
219 CHKERR set_dm_section(dm);
221 boost::shared_ptr<FEMethod>
null;
248 auto copy_dm_struture = [&](
auto simple_dm) {
259 dm = copy_dm_struture(
simple->getDM());
264 auto set_dm_section = [&](
auto dm) {
268 CHKERR DMSetSection(dm, section);
271 CHKERR set_dm_section(dm);
273 boost::shared_ptr<FEMethod>
null;
309 auto copy_dm_struture = [&](
auto simple_dm) {
320 dm = copy_dm_struture(
simple->getDM());
325 auto set_dm_section = [&](
auto dm) {
329 CHKERR DMSetSection(dm, section);
332 CHKERR set_dm_section(dm);
334 boost::shared_ptr<FEMethod>
null;
370 auto copy_dm_struture = [&](
auto simple_dm) {
381 dm = copy_dm_struture(
simple->getDM());
386 auto set_dm_section = [&](
auto dm) {
390 CHKERR DMSetSection(dm, section);
393 CHKERR set_dm_section(dm);
395 boost::shared_ptr<FEMethod>
null;
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set KSP right hand side evaluation function
SmartPetscObj< SNES > createSNES(SmartPetscObj< DM > dm=nullptr)
Create SNES (nonlinear) solver.
auto createTS(MPI_Comm comm)
virtual MPI_Comm & get_comm() const =0
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
auto createSNES(MPI_Comm comm)
auto getDMKspCtx(DM dm)
Get KSP context data structure used by DM.
MoFEMErrorCode loopFiniteElements(SmartPetscObj< DM > dm=nullptr)
Iterate finite elements.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
PetscObject getPetscObject(T obj)
boost::shared_ptr< FEMethod > feBoundaryRhs
Element to assemble RHS side by integrating boundary.
boost::shared_ptr< FEMethod > feDomainLhs
Element to assemble LHS side by integrating domain.
PipelineManager interface.
auto createKSP(MPI_Comm comm)
PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
Set KSP operators and push mofem finite element methods.
SmartPetscObj< TS > createTSIM(SmartPetscObj< DM > dm=nullptr)
Create TS (time) implicit solver.
Simple interface for fast problem set-up.
Deprecated interface functions.
boost::shared_ptr< FEMethod > feBoundaryLhs
Element to assemble LHS side by integrating boundary.
#define CHKERR
Inline error check.
PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS Jacobian evaluation function
implementation of Data Operators for Forces and Sources
Section manager is used to create indexes and sections.
auto createDM(MPI_Comm comm, const std::string dm_type_name)
Creates smart DM object.
PetscErrorCode DMMoFEMTSSetI2Jacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS Jacobian evaluation function
boost::shared_ptr< FEMethod > feDomainRhs
Element to assemble RHS side by integrating domain.
PetscErrorCode DMMoFEMTSSetI2Function(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS implicit function evaluation function
void simple(double P1[], double P2[], double P3[], double c[], const int N)
SmartPetscObj< KSP > createKSP(SmartPetscObj< DM > dm=nullptr)
Create KSP (linear) solver.
boost::shared_ptr< FEMethod > feBoundaryExplicitRhs
PipelineManager(const MoFEM::Core &core)
PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES Jacobian evaluation function
boost::shared_ptr< FEMethod > feSkeletonLhs
Element to assemble LHS side by integrating skeleton.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
base class for all interface classes
SmartPetscObj< TS > createTSEX(SmartPetscObj< DM > dm=nullptr)
Create TS (time) explit solver.
PetscErrorCode DMMoFEMTSSetRHSFunction(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS the right hand side function
PetscErrorCode DMMoFEMTSSetIFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set TS implicit function evaluation function
SmartPetscObj< TS > createTSIMEX(SmartPetscObj< DM > dm=nullptr)
Create TS (time) implicit-explicit solver.
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
boost::shared_ptr< FEMethod > feDomainExplicitRhs
Element to assemble explict Rhs for IMEX solver.
auto getDMSnesCtx(DM dm)
Get SNES context data structure used by DM.
keeps basic data about problem
SmartPetscObj< TS > createTS(const TSType type, SmartPetscObj< DM > dm=nullptr)
create TS (time) solver
PetscErrorCode DMMoFEMDuplicateDMCtx(DM dm, DM dm_duplicate)
Duplicate internal data struture.
boost::shared_ptr< FEMethod > feSkeletonRhs
Element to assemble RHS side by integrating skeleton.
SmartPetscObj< TS > createTSIM2(SmartPetscObj< DM > dm=nullptr)
Create TS (time) solver for second order equation in time.
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
auto getDMTsCtx(DM dm)
Get TS context data structure used by DM.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
boost::shared_ptr< FEMethod > feSkeletonExplicitRhs
PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES residual evaluation function