|
| v0.14.0
|
Go to the documentation of this file.
7 #ifndef __PCMGSETUP_VIA_APPROX_ORDERS_HPP__
8 #define __PCMGSETUP_VIA_APPROX_ORDERS_HPP__
35 boost::ptr_vector<PCMGSubMatrixCtx>
81 bool create_sub_matrix,
260 #endif //__PCMGSETUP_VIA_APPROX_ORDERS_HPP__
virtual ~PCMGSetUpViaApproxOrdersCtx()=default
MoFEMErrorCode DMMGViaApproxOrdersReplaceCoarseningIS(DM dm, IS *is_vec, int nb_elems, Mat A, int verb=0)
Replace coarsening IS in DM via approximation orders.
MoFEMErrorCode DMMGViaApproxOrdersClearCoarseningIS(DM)
Clear approximation orders.
boost::ptr_vector< PCMGSubMatrixCtx > shellMatrixCtxPtr
Shell sub-matrix context.
MoFEMErrorCode DMCoarsen_MGViaApproxOrders(DM dm, MPI_Comm comm, DM *dmc)
Coarsen DM.
std::vector< Mat > kspOperators
Get KSP operators.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Set data structures of MG pre-conditioner via approximation orders.
Structure for DM for multi-grid via approximation orders.
Implementation of DM context. You should not use it directly.
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.
int nbLevels
number of multi-grid levels
PCMGSetUpViaApproxOrdersCtx(DM dm, Mat a, bool shell_sub_a)
MoFEMErrorCode DMCreateInterpolation_MGViaApproxOrders(DM dm1, DM dm2, Mat *mat, Vec *vec)
Create interpolation matrix between data managers dm1 and dm2.
FTensor::Index< 'M', 3 > M
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, MoFEM::UnknownInterface **iface) const
DM dM
Distributed mesh manager.
virtual ~DMMGViaApproxOrdersCtx()
MoFEMErrorCode destroyCoarseningIS()
virtual MoFEMErrorCode createIsAtLevel(int kk, IS *is)
Set IS for levels.
MoFEMErrorCode DMCreateGlobalVector_MGViaApproxOrders(DM dm, Vec *g)
Create global vector for DMGViaApproxOrders.
MoFEMErrorCode DMCreate_MGViaApproxOrders(DM dm)
Create DM data structure for Multi-Grid via approximation orders.
MoFEMErrorCode DMMGViaApproxOrdersGetCoarseningISSize(DM dm, int *size)
Gets size of coarseningIS in internal data struture DMMGViaApproxOrdersCtx.
MoFEMErrorCode DMCreateMatrix_MGViaApproxOrders(DM dm, Mat *M)
Create matrix for Multi-Grid via approximation orders.
MoFEMErrorCode DMMGViaApproxOrdersPopBackCoarseningIS(DM)
Pop is form MG via approximation orders.
int orderAtLastLevel
set maximal evaluated order
MoFEMErrorCode PCMGSetUpViaApproxOrders(PC pc, PCMGSetUpViaApproxOrdersCtx *ctx, int verb=0)
Function build MG structure.
base class for all interface classes
MoFEMErrorCode DMRegister_MGViaApproxOrders(const char sname[])
Register DM for Multi-Grid via approximation orders.
virtual MoFEMErrorCode buildProlongationOperator(bool use_mat_a, int verb=0)
Set up data structures for MG.
std::vector< IS > coarseningIS
Coarsening IS.
MoFEMErrorCode DMMGViaApproxOrdersGetCtx(DM dm, DMMGViaApproxOrdersCtx **ctx)
int coarseOrder
approximation order of coarse level
virtual MoFEMErrorCode destroyIsAtLevel(int kk, IS *is)
Destroy IS if internally created.
const FTensor::Tensor2< T, Dim, Dim > Vec
PCMGSubMatrixCtx(Mat a, IS is)
virtual MoFEMErrorCode getOptions()
get options from line command
MoFEMErrorCode DMMGViaApproxOrdersSetAO(DM dm, AO ao)
Set DM ordering.
Mat A
Matrix at fine level.
virtual ~PCMGSubMatrixCtx()