38 #include "PetscVecTools.hpp" 39 #include "PCLDUFactorisation.hpp" 41 #include "Warnings.hpp" 60 PCLDUFactorisation::~PCLDUFactorisation()
65 std::cout <<
" -- LDU factorisation preconditioner profile information: " << std::endl;
66 std::cout <<
"\t mScatterTime: " <<
mPCContext.mScatterTime << std::endl;
67 std::cout <<
"\t mA1PreconditionerTime: " <<
mPCContext.mA1PreconditionerTime << std::endl;
68 std::cout <<
"\t mA2PreconditionerTime: " <<
mPCContext.mA2PreconditionerTime << std::endl;
69 std::cout <<
"\t mExtraLAOperations: " <<
mPCContext.mExtraLAOperations << std::endl;
70 std::cout <<
"\t mGatherTime: " <<
mPCContext.mGatherTime << std::endl;
96 Mat system_matrix, dummy;
97 #if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5) 98 KSPGetOperators(rKspObject, &system_matrix, &dummy);
101 KSPGetOperators(rKspObject, &system_matrix, &dummy, &flag);
105 MatGetSize(system_matrix, &num_rows, &num_columns);
107 PetscInt num_local_rows, num_local_columns;
108 MatGetLocalSize(system_matrix, &num_local_rows, &num_local_columns);
112 if ((num_rows%2 != 0) || (num_local_rows%2 != 0))
114 TERMINATE(
"Wrong matrix parallel layout detected in PCLDUFactorisation.");
118 unsigned subvector_num_rows = num_rows/2;
119 unsigned subvector_local_rows = num_local_rows/2;
143 assert(global_size == num_rows/2);
147 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &A11_local_rows);
148 ISCreateStride(PETSC_COMM_WORLD, global_size, 0, 2, &A11_columns);
150 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later 151 MatGetSubMatrix(system_matrix, A11_local_rows, A11_local_rows,
154 MatGetSubMatrix(system_matrix, A11_local_rows, A11_columns, PETSC_DECIDE,
167 assert(global_size == num_rows/2);
171 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low+1, 2, &A22_local_rows);
172 ISCreateStride(PETSC_COMM_WORLD, global_size, 1, 2, &A22_columns);
174 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later 175 MatGetSubMatrix(system_matrix, A22_local_rows, A22_local_rows,
178 MatGetSubMatrix(system_matrix, A22_local_rows, A22_columns, PETSC_DECIDE,
192 assert(global_size == num_rows/2);
196 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &B_local_rows);
198 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later 199 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low+1, 2, &B_columns);
200 MatGetSubMatrix(system_matrix, B_local_rows, B_columns,
203 ISCreateStride(PETSC_COMM_WORLD, global_size, 1, 2, &B_columns);
204 MatGetSubMatrix(system_matrix, B_local_rows, B_columns, PETSC_DECIDE,
233 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2 255 #if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5) 270 PetscPushErrorHandler(PetscIgnoreErrorHandler,
nullptr);
272 if (pc_set_error != 0)
274 WARNING(
"PETSc hypre preconditioning library is not installed");
277 PetscPopErrorHandler();
279 #if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR<=5) 315 #if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5) 328 PetscPushErrorHandler(PetscIgnoreErrorHandler,
nullptr);
331 PetscPopErrorHandler();
360 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later 361 PetscErrorCode PCLDUFactorisationApply(PC pc_object,
Vec x,
Vec y)
365 PCShellGetContext(pc_object, &pc_context);
367 PetscErrorCode PCLDUFactorisationApply(
void* pc_context,
Vec x,
Vec y)
374 assert(block_diag_context!=
nullptr);
410 double minus_one = -1.0;
411 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2 412 VecAYPX(&minus_one, block_diag_context->
x2_subvector, block_diag_context->
temp);
414 VecAYPX(block_diag_context->
temp, minus_one, block_diag_context->
x2_subvector);
448 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2 449 VecAYPX(&minus_one, block_diag_context->
z, block_diag_context->
y1_subvector);
451 VecAYPX(block_diag_context->
y1_subvector, minus_one, block_diag_context->
z);
void PCLDUFactorisationSetUp()
PCLDUFactorisationContext mPCContext
#define TERMINATE(message)
VecScatter A11_scatter_ctx
static double GetElapsedTime()
PCLDUFactorisation(KSP &rKspObject)
VecScatter A22_scatter_ctx
void PCLDUFactorisationCreate(KSP &rKspObject)