v0.14.0 |
Implementation of PETSc DM, managing interactions between mesh data structures and vectors and matrices. More...
Functions | |
PetscErrorCode | MoFEM::DMRegister_MoFEM (const char sname[]) |
Register MoFEM problem. More... | |
PetscErrorCode | MoFEM::DMMoFEMCreateMoFEM (DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set()) |
Must be called by user to set MoFEM data structures. More... | |
PetscErrorCode | MoFEM::DMMoFEMCreateSubDM (DM subdm, DM dm, const char problem_name[]) |
Must be called by user to set Sub DM MoFEM data structures. More... | |
PetscErrorCode | MoFEM::DMoFEMGetInterfacePtr (DM dm, MoFEM::Interface **m_field_ptr) |
Get pointer to MoFEM::Interface. More... | |
PetscErrorCode | MoFEM::DMMoFEMGetProblemPtr (DM dm, const MoFEM::Problem **problem_ptr) |
Get pointer to problem data structure. More... | |
PetscErrorCode | MoFEM::DMMoFEMSetSquareProblem (DM dm, PetscBool square_problem) |
set squared problem More... | |
PetscErrorCode | MoFEM::DMMoFEMGetSquareProblem (DM dm, PetscBool *square_problem) |
get squared problem More... | |
PetscErrorCode | MoFEM::DMMoFEMResolveSharedFiniteElements (DM dm, std::string fe_name) |
Resolve shared entities. More... | |
PetscErrorCode | MoFEM::DMMoFEMGetProblemFiniteElementLayout (DM dm, std::string fe_name, PetscLayout *layout) |
Get finite elements layout in the problem. More... | |
PetscErrorCode | MoFEM::DMMoFEMAddElement (DM dm, std::string fe_name) |
add element to dm More... | |
PetscErrorCode | MoFEM::DMMoFEMAddElement (DM dm, std::vector< std::string > fe_name) |
add element to dm More... | |
PetscErrorCode | MoFEM::DMMoFEMUnSetElement (DM dm, std::string fe_name) |
unset element from dm More... | |
PetscErrorCode | MoFEM::DMoFEMMeshToLocalVector (DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode) |
set local (or ghosted) vector values on mesh for partition only More... | |
PetscErrorCode | MoFEM::DMoFEMMeshToGlobalVector (DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode) |
set ghosted vector values on all existing mesh entities More... | |
PetscErrorCode | MoFEM::DMoFEMPreProcessFiniteElements (DM dm, MoFEM::FEMethod *method) |
execute finite element method for each element in dm (problem) More... | |
PetscErrorCode | MoFEM::DMoFEMPostProcessFiniteElements (DM dm, MoFEM::FEMethod *method) |
execute finite element method for each element in dm (problem) More... | |
PetscErrorCode | MoFEM::DMoFEMLoopFiniteElementsUpAndLowRank (DM dm, const char fe_name[], MoFEM::FEMethod *method, int low_rank, int up_rank, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr()) |
Executes FEMethod for finite elements in DM. More... | |
PetscErrorCode | MoFEM::DMoFEMLoopFiniteElementsUpAndLowRank (DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, int low_rank, int up_rank, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr()) |
Executes FEMethod for finite elements in DM. More... | |
PetscErrorCode | MoFEM::DMoFEMLoopFiniteElements (DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr()) |
Executes FEMethod for finite elements in DM. More... | |
PetscErrorCode | MoFEM::DMoFEMLoopFiniteElements (DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr()) |
Executes FEMethod for finite elements in DM. More... | |
PetscErrorCode | MoFEM::DMoFEMLoopDofs (DM dm, const char field_name[], MoFEM::DofMethod *method) |
execute method for dofs on field in problem More... | |
PetscErrorCode | MoFEM::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 More... | |
PetscErrorCode | MoFEM::DMMoFEMKSPSetComputeRHS (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 KSP right hand side evaluation function More... | |
PetscErrorCode | MoFEM::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. More... | |
PetscErrorCode | MoFEM::DMMoFEMKSPSetComputeOperators (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 KSP operators and push mofem finite element methods. More... | |
PetscErrorCode | MoFEM::DMMoFEMSNESSetFunction (DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only) |
set SNES residual evaluation function More... | |
PetscErrorCode | MoFEM::DMMoFEMSNESSetFunction (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 SNES residual evaluation function More... | |
PetscErrorCode | MoFEM::DMMoFEMSNESSetJacobian (DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only) |
set SNES Jacobian evaluation function More... | |
PetscErrorCode | MoFEM::DMMoFEMSNESSetJacobian (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 SNES Jacobian evaluation function More... | |
PetscErrorCode | MoFEM::DMMoFEMTSSetIFunction (DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only) |
set TS implicit function evaluation function More... | |
PetscErrorCode | MoFEM::DMMoFEMTSSetIFunction (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 More... | |
PetscErrorCode | MoFEM::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 More... | |
PetscErrorCode | MoFEM::DMMoFEMTSSetIJacobian (DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only) |
set TS Jacobian evaluation function More... | |
PetscErrorCode | MoFEM::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 More... | |
PetscErrorCode | MoFEM::DMMoFEMTSSetI2Function (DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only) |
set TS implicit function evaluation function More... | |
PetscErrorCode | MoFEM::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 More... | |
PetscErrorCode | MoFEM::DMMoFEMTSSetI2Jacobian (DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only) |
set TS Jacobian evaluation function More... | |
PetscErrorCode | MoFEM::DMMoFEMGetKspCtx (DM dm, MoFEM::KspCtx **ksp_ctx) |
get MoFEM::KspCtx data structure More... | |
PetscErrorCode | MoFEM::DMMoFEMGetKspCtx (DM dm, const boost::shared_ptr< MoFEM::KspCtx > &ksp_ctx) |
get MoFEM::KspCtx data structure More... | |
PetscErrorCode | MoFEM::DMMoFEMSetKspCtx (DM dm, boost::shared_ptr< MoFEM::KspCtx > ksp_ctx) |
set MoFEM::KspCtx data structure More... | |
PetscErrorCode | MoFEM::DMMoFEMGetSnesCtx (DM dm, MoFEM::SnesCtx **snes_ctx) |
get MoFEM::SnesCtx data structure More... | |
PetscErrorCode | MoFEM::DMMoFEMGetSnesCtx (DM dm, const boost::shared_ptr< MoFEM::SnesCtx > &snes_ctx) |
get MoFEM::SnesCtx data structure More... | |
PetscErrorCode | MoFEM::DMMoFEMSetSnesCtx (DM dm, boost::shared_ptr< MoFEM::SnesCtx > snes_ctx) |
Set MoFEM::SnesCtx data structure. More... | |
PetscErrorCode | MoFEM::DMMoFEMGetTsCtx (DM dm, MoFEM::TsCtx **ts_ctx) |
get MoFEM::TsCtx data structure More... | |
PetscErrorCode | MoFEM::DMMoFEMGetTsCtx (DM dm, const boost::shared_ptr< MoFEM::TsCtx > &ts_ctx) |
get MoFEM::TsCtx data structure More... | |
PetscErrorCode | MoFEM::DMMoFEMSetTsCtx (DM dm, boost::shared_ptr< MoFEM::TsCtx > ts_ctx) |
Set MoFEM::TsCtx data structure. More... | |
PetscErrorCode | MoFEM::DMMoFEMSetIsPartitioned (DM dm, PetscBool is_partitioned) |
PetscErrorCode | MoFEM::DMMoFEMGetIsPartitioned (DM dm, PetscBool *is_partitioned) |
PetscErrorCode | MoFEM::DMSetOperators_MoFEM (DM dm) |
Set operators for MoFEM dm. More... | |
PetscErrorCode | MoFEM::DMCreate_MoFEM (DM dm) |
Create dm data structure with MoFEM data structure. More... | |
PetscErrorCode | MoFEM::DMDestroy_MoFEM (DM dm) |
Destroys dm with MoFEM data structure. More... | |
PetscErrorCode | MoFEM::DMCreateGlobalVector_MoFEM (DM dm, Vec *g) |
DMShellSetCreateGlobalVector. More... | |
PetscErrorCode | MoFEM::DMCreateGlobalVector_MoFEM (DM dm, SmartPetscObj< Vec > &g_ptr) |
DMShellSetCreateGlobalVector. More... | |
PetscErrorCode | MoFEM::DMCreateLocalVector_MoFEM (DM dm, Vec *l) |
DMShellSetCreateLocalVector. More... | |
PetscErrorCode | MoFEM::DMCreateMatrix_MoFEM (DM dm, Mat *M) |
PetscErrorCode | MoFEM::DMCreateMatrix_MoFEM (DM dm, SmartPetscObj< Mat > &M) |
PetscErrorCode | MoFEM::DMSetFromOptions_MoFEM (DM dm) |
PetscErrorCode | MoFEM::DMSetUp_MoFEM (DM dm) |
PetscErrorCode | MoFEM::DMSubDMSetUp_MoFEM (DM subdm) |
PetscErrorCode | MoFEM::DMMoFEMAddSubFieldRow (DM dm, const char field_name[]) |
PetscErrorCode | MoFEM::DMMoFEMAddSubFieldRow (DM dm, const char field_name[], boost::shared_ptr< Range > r_ptr) |
Add field to sub dm problem on rows. More... | |
PetscErrorCode | MoFEM::DMMoFEMAddSubFieldCol (DM dm, const char field_name[]) |
PetscErrorCode | MoFEM::DMMoFEMGetIsSubDM (DM dm, PetscBool *is_sub_dm) |
PetscErrorCode | MoFEM::DMMoFEMAddRowCompositeProblem (DM dm, const char prb_name[]) |
Add problem to composite DM on row. More... | |
PetscErrorCode | MoFEM::DMMoFEMAddColCompositeProblem (DM dm, const char prb_name[]) |
Add problem to composite DM on col. More... | |
PetscErrorCode | MoFEM::DMMoFEMGetIsCompDM (DM dm, PetscBool *is_comp_dm) |
Get if this DM is composite DM. More... | |
PetscErrorCode | MoFEM::DMGlobalToLocalBegin_MoFEM (DM dm, Vec, InsertMode, Vec) |
PetscErrorCode | MoFEM::DMGlobalToLocalEnd_MoFEM (DM dm, Vec, InsertMode, Vec) |
PetscErrorCode | MoFEM::DMLocalToGlobalBegin_MoFEM (DM, Vec, InsertMode, Vec) |
PetscErrorCode | MoFEM::DMLocalToGlobalEnd_MoFEM (DM, Vec, InsertMode, Vec) |
PetscErrorCode | MoFEM::DMCreateFieldIS_MoFEM (DM dm, PetscInt *numFields, char ***fieldNames, IS **fields) |
PetscErrorCode | MoFEM::DMMoFEMGetFieldIS (DM dm, RowColData rc, const char field_name[], IS *is) |
get field is in the problem More... | |
auto | MoFEM::createDMMatrix (DM dm) |
Get smart matrix from DM. More... | |
auto | MoFEM::createDMVector (DM dm) |
Get smart vector from DM. More... | |
MoFEMErrorCode | DMMGViaApproxOrdersPushBackCoarseningIS (DM, IS is, Mat A, Mat *subA, bool create_sub_matrix, bool shell_sub_a) |
Push back coarsening level to MG via approximation orders. More... | |
MoFEMErrorCode | DMMGViaApproxOrdersPopBackCoarseningIS (DM) |
Pop is form MG via approximation orders. More... | |
MoFEMErrorCode | DMMGViaApproxOrdersClearCoarseningIS (DM) |
Clear approximation orders. More... | |
MoFEMErrorCode | DMRegister_MGViaApproxOrders (const char sname[]) |
Register DM for Multi-Grid via approximation orders. More... | |
MoFEMErrorCode | DMCreateMatrix_MGViaApproxOrders (DM dm, Mat *M) |
Create matrix for Multi-Grid via approximation orders. More... | |
MoFEMErrorCode | DMCoarsen_MGViaApproxOrders (DM dm, MPI_Comm comm, DM *dmc) |
Coarsen DM. More... | |
Implementation of PETSc DM, managing interactions between mesh data structures and vectors and matrices.
DM objects are used to manage communication between the algebraic structures in PETSc (Vec and Mat) and mesh data structures in PDE-based (or other) simulations.
DM is abstract interface, here is it particular implementation for MoFEM code.
|
inline |
Get smart matrix from DM.
Definition at line 1003 of file DMMoFEM.hpp.
|
inline |
Get smart vector from DM.
Definition at line 1018 of file DMMoFEM.hpp.
MoFEMErrorCode DMCoarsen_MGViaApproxOrders | ( | DM | dm, |
MPI_Comm | comm, | ||
DM * | dmc | ||
) |
Coarsen DM.
Not used directly by user
dm | Distributed mesh data structure |
comm | Communicator |
dmc | Coarse distributed mesh data structure |
Definition at line 462 of file PCMGSetUpViaApproxOrders.cpp.
PetscErrorCode MoFEM::DMCreate_MoFEM | ( | DM | dm | ) |
PetscErrorCode MoFEM::DMCreateFieldIS_MoFEM | ( | DM | dm, |
PetscInt * | numFields, | ||
char *** | fieldNames, | ||
IS ** | fields | ||
) |
Creates a set of IS objects with the global indices of dofs for each field
dm | The number of fields (or NULL if not requested) |
Output:
numFields | The number of fields (or NULL if not requested) |
fieldNames | The name for each field (or NULL if not requested) |
fields | The global indices for each field (or NULL if not requested) |
Definition at line 1451 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMCreateGlobalVector_MoFEM | ( | DM | dm, |
SmartPetscObj< Vec > & | g_ptr | ||
) |
DMShellSetCreateGlobalVector.
sets the routine to create a global vector associated with the shell DM
Definition at line 1181 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMCreateGlobalVector_MoFEM | ( | DM | dm, |
Vec * | g | ||
) |
DMShellSetCreateGlobalVector.
sets the routine to create a global vector associated with the shell DM
Definition at line 1171 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMCreateLocalVector_MoFEM | ( | DM | dm, |
Vec * | l | ||
) |
DMShellSetCreateLocalVector.
sets the routine to create a local vector associated with the shell DM
Definition at line 1191 of file DMMoFEM.cpp.
MoFEMErrorCode DMCreateMatrix_MGViaApproxOrders | ( | DM | dm, |
Mat * | M | ||
) |
Create matrix for Multi-Grid via approximation orders.
Not used directly by user
dm | Distributed mesh data structure |
M | Matrix |
Definition at line 430 of file PCMGSetUpViaApproxOrders.cpp.
PetscErrorCode MoFEM::DMCreateMatrix_MoFEM | ( | DM | dm, |
Mat * | M | ||
) |
DMShellSetCreateMatrix
sets the routine to create a matrix associated with the shell DM
Definition at line 1201 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMCreateMatrix_MoFEM | ( | DM | dm, |
SmartPetscObj< Mat > & | M | ||
) |
DMShellSetCreateMatrix
sets the routine to create a matrix associated with the shell DM
Definition at line 1229 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMDestroy_MoFEM | ( | DM | dm | ) |
Destroys dm with MoFEM data structure.
destroy the MoFEM structure
Definition at line 83 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMGlobalToLocalBegin_MoFEM | ( | DM | dm, |
Vec | g, | ||
InsertMode | mode, | ||
Vec | l | ||
) |
DMShellSetGlobalToLocal
the routine that begins the global to local scatter
Definition at line 1379 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMGlobalToLocalEnd_MoFEM | ( | DM | dm, |
Vec | g, | ||
InsertMode | mode, | ||
Vec | l | ||
) |
DMShellSetGlobalToLocal
the routine that begins the global to local scatter
Definition at line 1387 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMLocalToGlobalBegin_MoFEM | ( | DM | dm, |
Vec | l, | ||
InsertMode | mode, | ||
Vec | g | ||
) |
DMShellSetLocalToGlobal
the routine that begins the local to global scatter
DMShellSetLocalToGlobal
the routine that ends the local to global scatter
Definition at line 1416 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMLocalToGlobalEnd_MoFEM | ( | DM | , |
Vec | l, | ||
InsertMode | mode, | ||
Vec | g | ||
) |
DMShellSetLocalToGlobal
the routine that ends the local to global scatter
Definition at line 1445 of file DMMoFEM.cpp.
MoFEMErrorCode DMMGViaApproxOrdersClearCoarseningIS | ( | DM | ) |
Clear approximation orders.
DM | dm |
Definition at line 276 of file PCMGSetUpViaApproxOrders.cpp.
MoFEMErrorCode DMMGViaApproxOrdersPopBackCoarseningIS | ( | DM | ) |
Pop is form MG via approximation orders.
DM | dm |
is | pop back IS |
Definition at line 260 of file PCMGSetUpViaApproxOrders.cpp.
MoFEMErrorCode DMMGViaApproxOrdersPushBackCoarseningIS | ( | DM | , |
IS | is, | ||
Mat | A, | ||
Mat * | subA, | ||
bool | create_sub_matrix, | ||
bool | shell_sub_a | ||
) |
Push back coarsening level to MG via approximation orders.
DM | discrete manager |
is | Push back IS used for coarsening |
A | Get sub-matrix of A using is (that sets operators for coarsening levels) |
subA | Returning pointer to created sub matrix |
subA | If true create sub matrix, otherwise in subA has to be valid pointer to subA |
Definition at line 206 of file PCMGSetUpViaApproxOrders.cpp.
PetscErrorCode MoFEM::DMMoFEMAddColCompositeProblem | ( | DM | dm, |
const char | prb_name[] | ||
) |
Add problem to composite DM on col.
This create block on col with DOFs from problem of given name
dm | the DM object |
prb_name | add problem name |
Definition at line 389 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMAddElement | ( | DM | dm, |
std::string | fe_name | ||
) |
add element to dm
Definition at line 501 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMAddElement | ( | DM | dm, |
std::vector< std::string > | fe_name | ||
) |
add element to dm
Definition at line 510 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMAddRowCompositeProblem | ( | DM | dm, |
const char | prb_name[] | ||
) |
Add problem to composite DM on row.
This create block on row with DOFs from problem of given name
dm | the DM object |
prb_name | add problem name |
Definition at line 371 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMAddSubFieldCol | ( | DM | dm, |
const char | field_name[] | ||
) |
Add field to sub dm problem on columns
Definition at line 284 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMAddSubFieldRow | ( | DM | dm, |
const char | field_name[] | ||
) |
Add field to sub dm problem on rows
Definition at line 242 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMAddSubFieldRow | ( | DM | dm, |
const char | field_name[], | ||
boost::shared_ptr< Range > | r_ptr | ||
) |
Add field to sub dm problem on rows.
dm | |
field_name | |
m |
Definition at line 262 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMCreateMoFEM | ( | DM | dm, |
MoFEM::Interface * | m_field_ptr, | ||
const char | problem_name[], | ||
const MoFEM::BitRefLevel | bit_level, | ||
const MoFEM::BitRefLevel | bit_mask = MoFEM::BitRefLevel().set() |
||
) |
Must be called by user to set MoFEM data structures.
Definition at line 118 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMCreateSubDM | ( | DM | subdm, |
DM | dm, | ||
const char | problem_name[] | ||
) |
Must be called by user to set Sub DM MoFEM data structures.
Definition at line 219 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMGetFieldIS | ( | DM | dm, |
RowColData | rc, | ||
const char | field_name[], | ||
IS * | is | ||
) |
get field is in the problem
dm | the DM objects |
rc | ROW or COL (no difference is problem is squared) |
field_name | name of the field |
is | returned the IS object |
Definition at line 1496 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMGetIsCompDM | ( | DM | dm, |
PetscBool * | is_comp_dm | ||
) |
Get if this DM is composite DM.
dm | the DM object |
is_comp_dm | return true if composite problem here |
Definition at line 409 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMGetIsPartitioned | ( | DM | dm, |
PetscBool * | is_partitioned | ||
) |
PetscErrorCode MoFEM::DMMoFEMGetIsSubDM | ( | DM | dm, |
PetscBool * | is_sub_dm | ||
) |
Return true if this DM is sub problem
dm | the DM object |
is_subproblem | true if subproblem |
Definition at line 326 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMGetKspCtx | ( | DM | dm, |
const boost::shared_ptr< MoFEM::KspCtx > & | ksp_ctx | ||
) |
PetscErrorCode MoFEM::DMMoFEMGetKspCtx | ( | DM | dm, |
MoFEM::KspCtx ** | ksp_ctx | ||
) |
PetscErrorCode MoFEM::DMMoFEMGetProblemFiniteElementLayout | ( | DM | dm, |
std::string | fe_name, | ||
PetscLayout * | layout | ||
) |
Get finite elements layout in the problem.
In layout is stored information how many elements is on each processor, for more information look int petsc documentation http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/PetscLayoutCreate.html#PetscLayoutCreate
dm | discrete manager for this problem |
fe_name | finite element name |
layout | pointer to layout, for created layout user takes responsibility for destroying it. |
Definition at line 477 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMGetProblemPtr | ( | DM | dm, |
const MoFEM::Problem ** | problem_ptr | ||
) |
Get pointer to problem data structure.
Definition at line 430 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMGetSnesCtx | ( | DM | dm, |
const boost::shared_ptr< MoFEM::SnesCtx > & | snes_ctx | ||
) |
PetscErrorCode MoFEM::DMMoFEMGetSnesCtx | ( | DM | dm, |
MoFEM::SnesCtx ** | snes_ctx | ||
) |
get MoFEM::SnesCtx data structure
Definition at line 1098 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMGetSquareProblem | ( | DM | dm, |
PetscBool * | square_problem | ||
) |
get squared problem
It if true is assumed that matrix has the same indexing on rows and columns. This reduces interprocessor communication.
Definition at line 491 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMGetTsCtx | ( | DM | dm, |
const boost::shared_ptr< MoFEM::TsCtx > & | ts_ctx | ||
) |
PetscErrorCode MoFEM::DMMoFEMGetTsCtx | ( | DM | dm, |
MoFEM::TsCtx ** | ts_ctx | ||
) |
get MoFEM::TsCtx data structure
Definition at line 1146 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::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.
dm | DM |
fe_name | finite element name |
method | method on the element (executed for each element in the problem which given name) |
pre_only | method for pre-process before element method |
post_only | method for post-process after element method |
Definition at line 682 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMKSPSetComputeOperators | ( | 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 KSP operators and push mofem finite element methods.
dm | DM |
fe_name | finite element name |
method | method on the element (executed for each element in the problem which given name) |
pre_only | method for pre-process before element method |
post_only | method for post-process after element method |
Definition at line 693 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::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
Definition at line 641 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMKSPSetComputeRHS | ( | 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 | ||
) |
PetscErrorCode MoFEM::DMMoFEMResolveSharedFiniteElements | ( | DM | dm, |
std::string | fe_name | ||
) |
Resolve shared entities.
dm | dm |
fe_name | finite element for which shared entities are resolved |
This allows for tag reduction or tag exchange, f.e.
Definition at line 468 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMSetIsPartitioned | ( | DM | dm, |
PetscBool | is_partitioned | ||
) |
sets if read mesh is partitioned
get if read mesh is partitioned
Definition at line 1127 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMSetKspCtx | ( | DM | dm, |
boost::shared_ptr< MoFEM::KspCtx > | ksp_ctx | ||
) |
PetscErrorCode MoFEM::DMMoFEMSetSnesCtx | ( | DM | dm, |
boost::shared_ptr< MoFEM::SnesCtx > | snes_ctx | ||
) |
PetscErrorCode MoFEM::DMMoFEMSetSquareProblem | ( | DM | dm, |
PetscBool | square_problem | ||
) |
set squared problem
It if true is assumed that matrix has the same indexing on rows and columns. This reduces interprocessor communication.
Definition at line 460 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMSetTsCtx | ( | DM | dm, |
boost::shared_ptr< MoFEM::TsCtx > | ts_ctx | ||
) |
Set MoFEM::TsCtx data structure.
It take over pointer, do not delete it, DM will destroy pointer when is destroyed.
Definition at line 1163 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMSNESSetFunction | ( | DM | dm, |
const char | fe_name[], | ||
MoFEM::FEMethod * | method, | ||
MoFEM::BasicMethod * | pre_only, | ||
MoFEM::BasicMethod * | post_only | ||
) |
set SNES residual evaluation function
Definition at line 722 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMSNESSetFunction | ( | 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 | ||
) |
PetscErrorCode MoFEM::DMMoFEMSNESSetJacobian | ( | DM | dm, |
const char | fe_name[], | ||
MoFEM::FEMethod * | method, | ||
MoFEM::BasicMethod * | pre_only, | ||
MoFEM::BasicMethod * | post_only | ||
) |
set SNES Jacobian evaluation function
Definition at line 763 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMSNESSetJacobian | ( | 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 | ||
) |
PetscErrorCode MoFEM::DMMoFEMTSSetI2Function | ( | DM | dm, |
const char | fe_name[], | ||
MoFEM::FEMethod * | method, | ||
MoFEM::BasicMethod * | pre_only, | ||
MoFEM::BasicMethod * | post_only | ||
) |
PetscErrorCode MoFEM::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
Definition at line 979 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMTSSetI2Jacobian | ( | DM | dm, |
const char | fe_name[], | ||
MoFEM::FEMethod * | method, | ||
MoFEM::BasicMethod * | pre_only, | ||
MoFEM::BasicMethod * | post_only | ||
) |
PetscErrorCode MoFEM::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
Definition at line 1021 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMTSSetIFunction | ( | DM | dm, |
const char | fe_name[], | ||
MoFEM::FEMethod * | method, | ||
MoFEM::BasicMethod * | pre_only, | ||
MoFEM::BasicMethod * | post_only | ||
) |
set TS implicit function evaluation function
Definition at line 804 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMTSSetIFunction | ( | 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 | ||
) |
PetscErrorCode MoFEM::DMMoFEMTSSetIJacobian | ( | DM | dm, |
const char | fe_name[], | ||
MoFEM::FEMethod * | method, | ||
MoFEM::BasicMethod * | pre_only, | ||
MoFEM::BasicMethod * | post_only | ||
) |
PetscErrorCode MoFEM::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
Definition at line 857 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMMoFEMUnSetElement | ( | DM | dm, |
std::string | fe_name | ||
) |
PetscErrorCode MoFEM::DMoFEMGetInterfacePtr | ( | DM | dm, |
MoFEM::Interface ** | m_field_ptr | ||
) |
Get pointer to MoFEM::Interface.
dm | Distributed mesh manager |
m_field_ptr | Pointer to pointer of field interface |
Definition at line 418 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMoFEMLoopDofs | ( | DM | dm, |
const char | field_name[], | ||
MoFEM::DofMethod * | method | ||
) |
PetscErrorCode MoFEM::DMoFEMLoopFiniteElements | ( | DM | dm, |
const char | fe_name[], | ||
MoFEM::FEMethod * | method, | ||
CacheTupleWeakPtr | cache_ptr = CacheTupleSharedPtr() |
||
) |
Executes FEMethod for finite elements in DM.
dm | MoFEM discrete manager |
fe_name | name of element |
method | pointer to MOFEM::FEMethod |
Definition at line 590 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMoFEMLoopFiniteElements | ( | DM | dm, |
const std::string | fe_name, | ||
boost::shared_ptr< MoFEM::FEMethod > | method, | ||
CacheTupleWeakPtr | cache_ptr = CacheTupleSharedPtr() |
||
) |
Executes FEMethod for finite elements in DM.
dm | MoFEM discrete manager |
fe_name | name of element |
method | pointer to MOFEM::FEMethod |
Definition at line 603 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMoFEMLoopFiniteElementsUpAndLowRank | ( | DM | dm, |
const char | fe_name[], | ||
MoFEM::FEMethod * | method, | ||
int | low_rank, | ||
int | up_rank, | ||
CacheTupleWeakPtr | cache_ptr = CacheTupleSharedPtr() |
||
) |
Executes FEMethod for finite elements in DM.
dm | MoFEM discrete manager |
fe_name | name of finite element |
method | pointer to MoFEM::FEMethod |
low_rank | lowest rank of processor |
up_rank | upper run of processor |
Definition at line 571 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMoFEMLoopFiniteElementsUpAndLowRank | ( | DM | dm, |
const std::string | fe_name, | ||
boost::shared_ptr< MoFEM::FEMethod > | method, | ||
int | low_rank, | ||
int | up_rank, | ||
CacheTupleWeakPtr | cache_ptr = CacheTupleSharedPtr() |
||
) |
Executes FEMethod for finite elements in DM.
dm | MoFEM discrete manager |
fe_name | name of finite element |
method | pointer to MoFEM::FEMethod |
low_rank | lowest rank of processor |
up_rank | upper run of processor |
Definition at line 583 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMoFEMMeshToGlobalVector | ( | DM | dm, |
Vec | g, | ||
InsertMode | mode, | ||
ScatterMode | scatter_mode | ||
) |
set ghosted vector values on all existing mesh entities
g | vector |
mode | see petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES) |
scatter_mode | see petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE) |
SCATTER_REVERSE set data to field entities from V vector.
SCATTER_FORWARD set vector V from data field entities
Definition at line 539 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMoFEMMeshToLocalVector | ( | DM | dm, |
Vec | l, | ||
InsertMode | mode, | ||
ScatterMode | scatter_mode | ||
) |
set local (or ghosted) vector values on mesh for partition only
l | vector |
mode | see petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES) |
scatter_mode | see petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE) |
SCATTER_REVERSE set data to field entities from V vector.
SCATTER_FORWARD set vector V from data field entities
Definition at line 527 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMoFEMPostProcessFiniteElements | ( | DM | dm, |
MoFEM::FEMethod * | method | ||
) |
execute finite element method for each element in dm (problem)
Definition at line 560 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMoFEMPreProcessFiniteElements | ( | DM | dm, |
MoFEM::FEMethod * | method | ||
) |
execute finite element method for each element in dm (problem)
Definition at line 550 of file DMMoFEM.cpp.
MoFEMErrorCode DMRegister_MGViaApproxOrders | ( | const char | sname[] | ) |
Register DM for Multi-Grid via approximation orders.
sname | problem/dm registered name |
Definition at line 376 of file PCMGSetUpViaApproxOrders.cpp.
PetscErrorCode MoFEM::DMRegister_MoFEM | ( | const char | sname[] | ) |
Register MoFEM problem.
Definition at line 47 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMSetFromOptions_MoFEM | ( | DM | dm | ) |
PetscErrorCode MoFEM::DMSetOperators_MoFEM | ( | DM | dm | ) |
Set operators for MoFEM dm.
dm |
Definition at line 53 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMSetUp_MoFEM | ( | DM | dm | ) |
sets up the MoFEM structures inside a DM object
Definition at line 1285 of file DMMoFEM.cpp.
PetscErrorCode MoFEM::DMSubDMSetUp_MoFEM | ( | DM | subdm | ) |
Sets up the MoFEM structures inside a DM object for sub dm
Definition at line 1332 of file DMMoFEM.cpp.