00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036 #include <iostream>
00037
00038 #include "PetscVecTools.hpp"
00039 #include "PCLDUFactorisation.hpp"
00040 #include "Exception.hpp"
00041 #include "Warnings.hpp"
00042
00043 PCLDUFactorisation::PCLDUFactorisation(KSP& rKspObject)
00044 {
00045 #ifdef TRACE_KSP
00046 mPCContext.mScatterTime = 0.0;
00047 mPCContext.mA1PreconditionerTime = 0.0;;
00048 mPCContext.mA2PreconditionerTime = 0.0;;
00049 mPCContext.mGatherTime = 0.0;;
00050 #endif
00051
00052 PCLDUFactorisationCreate(rKspObject);
00053 PCLDUFactorisationSetUp();
00054 }
00055
00056 PCLDUFactorisation::~PCLDUFactorisation()
00057 {
00058 #ifdef TRACE_KSP
00059 if (PetscTools::AmMaster())
00060 {
00061 std::cout << " -- LDU factorisation preconditioner profile information: " << std::endl;
00062 std::cout << "\t mScatterTime: " << mPCContext.mScatterTime << std::endl;
00063 std::cout << "\t mA1PreconditionerTime: " << mPCContext.mA1PreconditionerTime << std::endl;
00064 std::cout << "\t mA2PreconditionerTime: " << mPCContext.mA2PreconditionerTime << std::endl;
00065 std::cout << "\t mExtraLAOperations: " << mPCContext.mExtraLAOperations << std::endl;
00066 std::cout << "\t mGatherTime: " << mPCContext.mGatherTime << std::endl;
00067 }
00068 #endif
00069
00070 PetscTools::Destroy(mPCContext.A11_matrix_subblock);
00071 PetscTools::Destroy(mPCContext.A22_matrix_subblock);
00072 PetscTools::Destroy(mPCContext.B_matrix_subblock);
00073
00074 PCDestroy(PETSC_DESTROY_PARAM(mPCContext.PC_amg_A11));
00075 PCDestroy(PETSC_DESTROY_PARAM(mPCContext.PC_amg_A22));
00076
00077 PetscTools::Destroy(mPCContext.x1_subvector);
00078 PetscTools::Destroy(mPCContext.y1_subvector);
00079 PetscTools::Destroy(mPCContext.x2_subvector);
00080 PetscTools::Destroy(mPCContext.y2_subvector);
00081 PetscTools::Destroy(mPCContext.z);
00082 PetscTools::Destroy(mPCContext.temp);
00083
00084 VecScatterDestroy(PETSC_DESTROY_PARAM(mPCContext.A11_scatter_ctx));
00085 VecScatterDestroy(PETSC_DESTROY_PARAM(mPCContext.A22_scatter_ctx));
00086 }
00087
00088 void PCLDUFactorisation::PCLDUFactorisationCreate(KSP& rKspObject)
00089 {
00090 KSPGetPC(rKspObject, &mPetscPCObject);
00091
00092 Mat system_matrix, dummy;
00093 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
00094 KSPGetOperators(rKspObject, &system_matrix, &dummy);
00095 #else
00096 MatStructure flag;
00097 KSPGetOperators(rKspObject, &system_matrix, &dummy, &flag);
00098 #endif
00099
00100 PetscInt num_rows, num_columns;
00101 MatGetSize(system_matrix, &num_rows, &num_columns);
00102
00103 PetscInt num_local_rows, num_local_columns;
00104 MatGetLocalSize(system_matrix, &num_local_rows, &num_local_columns);
00105
00106
00107
00108 if ((num_rows%2 != 0) || (num_local_rows%2 != 0))
00109 {
00110 TERMINATE("Wrong matrix parallel layout detected in PCLDUFactorisation.");
00111 }
00112
00113
00114 unsigned subvector_num_rows = num_rows/2;
00115 unsigned subvector_local_rows = num_local_rows/2;
00116 mPCContext.x1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00117 mPCContext.x2_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00118 mPCContext.y1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00119 mPCContext.y2_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00120 mPCContext.z = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00121 mPCContext.temp = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00122
00123
00124 {
00125
00126 Vec dummy_vec = PetscTools::CreateVec(num_rows, num_local_rows);
00127
00128 PetscVecTools::SetupInterleavedVectorScatterGather(dummy_vec, mPCContext.A11_scatter_ctx, mPCContext.A22_scatter_ctx);
00129
00130 PetscTools::Destroy(dummy_vec);
00131 }
00132
00133
00134 {
00135
00136 PetscInt low, high, global_size;
00137 VecGetOwnershipRange(mPCContext.x1_subvector, &low, &high);
00138 VecGetSize(mPCContext.x1_subvector, &global_size);
00139 assert(global_size == num_rows/2);
00140
00141 IS A11_local_rows;
00142 IS A11_columns;
00143 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &A11_local_rows);
00144 ISCreateStride(PETSC_COMM_WORLD, global_size, 0, 2, &A11_columns);
00145
00146 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00147 MatGetSubMatrix(system_matrix, A11_local_rows, A11_local_rows,
00148 MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00149 #else
00150 MatGetSubMatrix(system_matrix, A11_local_rows, A11_columns, PETSC_DECIDE,
00151 MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00152 #endif
00153 ISDestroy(PETSC_DESTROY_PARAM(A11_local_rows));
00154 ISDestroy(PETSC_DESTROY_PARAM(A11_columns));
00155 }
00156
00157
00158 {
00159
00160 PetscInt low, high, global_size;
00161 VecGetOwnershipRange(mPCContext.x2_subvector, &low, &high);
00162 VecGetSize(mPCContext.x2_subvector, &global_size);
00163 assert(global_size == num_rows/2);
00164
00165 IS A22_local_rows;
00166 IS A22_columns;
00167 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low+1, 2, &A22_local_rows);
00168 ISCreateStride(PETSC_COMM_WORLD, global_size, 1, 2, &A22_columns);
00169
00170 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00171 MatGetSubMatrix(system_matrix, A22_local_rows, A22_local_rows,
00172 MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00173 #else
00174 MatGetSubMatrix(system_matrix, A22_local_rows, A22_columns, PETSC_DECIDE,
00175 MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00176 #endif
00177
00178 ISDestroy(PETSC_DESTROY_PARAM(A22_local_rows));
00179 ISDestroy(PETSC_DESTROY_PARAM(A22_columns));
00180 }
00181
00182
00183 {
00184
00185 PetscInt low, high, global_size;
00186 VecGetOwnershipRange(mPCContext.x1_subvector, &low, &high);
00187 VecGetSize(mPCContext.x1_subvector, &global_size);
00188 assert(global_size == num_rows/2);
00189
00190 IS B_local_rows;
00191 IS B_columns;
00192 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &B_local_rows);
00193
00194 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00195 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low+1, 2, &B_columns);
00196 MatGetSubMatrix(system_matrix, B_local_rows, B_columns,
00197 MAT_INITIAL_MATRIX, &mPCContext.B_matrix_subblock);
00198 #else
00199 ISCreateStride(PETSC_COMM_WORLD, global_size, 1, 2, &B_columns);
00200 MatGetSubMatrix(system_matrix, B_local_rows, B_columns, PETSC_DECIDE,
00201 MAT_INITIAL_MATRIX, &mPCContext.B_matrix_subblock);
00202 #endif
00203
00204 ISDestroy(PETSC_DESTROY_PARAM(B_local_rows));
00205 ISDestroy(PETSC_DESTROY_PARAM(B_columns));
00206 }
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219
00220
00221
00222
00223
00224
00225
00226
00227
00228 PCSetType(mPetscPCObject, PCSHELL);
00229 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00230
00231 PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply, (void*) &mPCContext);
00232 #else
00233
00234 PCShellSetContext(mPetscPCObject, &mPCContext);
00235
00236 PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply);
00237 #endif
00238 }
00239
00240 void PCLDUFactorisation::PCLDUFactorisationSetUp()
00241 {
00242
00243
00244
00245
00246
00247
00248
00249
00250 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A11));
00251 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
00252
00253 PCSetReusePreconditioner(mPCContext.PC_amg_A11, PETSC_TRUE);
00254 PCSetOperators(mPCContext.PC_amg_A11, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock);
00255 #else
00256 PCSetOperators(mPCContext.PC_amg_A11, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock, SAME_PRECONDITIONER);
00257 #endif
00258
00259
00260
00262
00263
00264
00265
00266 PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL);
00267 PetscErrorCode pc_set_error = PCSetType(mPCContext.PC_amg_A11, PCHYPRE);
00268 if (pc_set_error != 0)
00269 {
00270 WARNING("PETSc hypre preconditioning library is not installed");
00271 }
00272
00273 PetscPopErrorHandler();
00274
00275
00276 PetscOptionsSetValue("-pc_hypre_type", "euclid");
00277 PetscOptionsSetValue("-pc_hypre_euclid_levels", "0");
00278
00279
00280
00281
00282
00283
00284
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00302
00303 PCSetFromOptions(mPCContext.PC_amg_A11);
00304 PCSetUp(mPCContext.PC_amg_A11);
00305
00306
00307
00308
00309 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A22));
00310 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
00311
00312 PCSetReusePreconditioner(mPCContext.PC_amg_A22, PETSC_TRUE);
00313 PCSetOperators(mPCContext.PC_amg_A22, mPCContext.A22_matrix_subblock, mPCContext.A22_matrix_subblock);
00314 #else
00315 PCSetOperators(mPCContext.PC_amg_A22, mPCContext.A22_matrix_subblock, mPCContext.A22_matrix_subblock, SAME_PRECONDITIONER);
00316 #endif
00317
00318
00319
00321
00322
00323 PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL);
00324 PCSetType(mPCContext.PC_amg_A22, PCHYPRE);
00325
00326 PetscPopErrorHandler();
00327
00328 PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00329 PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00330 PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00331 PetscOptionsSetValue("-pc_hypre_boomeramg_coarsen_type", "HMIS");
00332
00333
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344
00345
00346
00347
00348
00349
00351
00352 PCSetFromOptions(mPCContext.PC_amg_A22);
00353 PCSetUp(mPCContext.PC_amg_A22);
00354 }
00355 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00356 PetscErrorCode PCLDUFactorisationApply(PC pc_object, Vec x, Vec y)
00357 {
00358 void* pc_context;
00359
00360 PCShellGetContext(pc_object, &pc_context);
00361 #else
00362 PetscErrorCode PCLDUFactorisationApply(void* pc_context, Vec x, Vec y)
00363 {
00364 #endif
00366
00367
00368 PCLDUFactorisation::PCLDUFactorisationContext* block_diag_context = (PCLDUFactorisation::PCLDUFactorisationContext*) pc_context;
00369 assert(block_diag_context!=NULL);
00370
00371
00372
00373
00374 #ifdef TRACE_KSP
00375 Timer::Reset();
00376 #endif
00377
00378 PetscVecTools::DoInterleavedVecScatter(x, block_diag_context->A11_scatter_ctx, block_diag_context->x1_subvector, block_diag_context->A22_scatter_ctx, block_diag_context->x2_subvector);
00379
00380 #ifdef TRACE_KSP
00381 block_diag_context->mScatterTime += Timer::GetElapsedTime();
00382 #endif
00383
00384
00385
00386
00387
00388
00389
00390
00391 #ifdef TRACE_KSP
00392 Timer::Reset();
00393 #endif
00394
00395 PCApply(block_diag_context->PC_amg_A11, block_diag_context->x1_subvector, block_diag_context->z);
00396 #ifdef TRACE_KSP
00397 block_diag_context->mA1PreconditionerTime += Timer::GetElapsedTime();
00398 #endif
00399
00400
00401 #ifdef TRACE_KSP
00402 Timer::Reset();
00403 #endif
00404 MatMult(block_diag_context->B_matrix_subblock,block_diag_context->z,block_diag_context->temp);
00405 double minus_one = -1.0;
00406 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00407 VecAYPX(&minus_one, block_diag_context->x2_subvector, block_diag_context->temp);
00408 #else
00409 VecAYPX(block_diag_context->temp, minus_one, block_diag_context->x2_subvector);
00410 #endif
00411 #ifdef TRACE_KSP
00412 block_diag_context->mExtraLAOperations += Timer::GetElapsedTime();
00413 #endif
00414
00415
00416 #ifdef TRACE_KSP
00417 Timer::Reset();
00418 #endif
00419 PCApply(block_diag_context->PC_amg_A22, block_diag_context->temp, block_diag_context->y2_subvector);
00420 #ifdef TRACE_KSP
00421 block_diag_context->mA2PreconditionerTime += Timer::GetElapsedTime();
00422 #endif
00423
00424
00425 #ifdef TRACE_KSP
00426 Timer::Reset();
00427 #endif
00428 MatMult(block_diag_context->B_matrix_subblock,block_diag_context->y2_subvector,block_diag_context->temp);
00429 #ifdef TRACE_KSP
00430 block_diag_context->mExtraLAOperations += Timer::GetElapsedTime();
00431 #endif
00432 #ifdef TRACE_KSP
00433 Timer::Reset();
00434 #endif
00435 PCApply(block_diag_context->PC_amg_A11, block_diag_context->temp, block_diag_context->y1_subvector);
00436 #ifdef TRACE_KSP
00437 block_diag_context->mA1PreconditionerTime += Timer::GetElapsedTime();
00438 #endif
00439
00440 #ifdef TRACE_KSP
00441 Timer::Reset();
00442 #endif
00443 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00444 VecAYPX(&minus_one, block_diag_context->z, block_diag_context->y1_subvector);
00445 #else
00446 VecAYPX(block_diag_context->y1_subvector, minus_one, block_diag_context->z);
00447 #endif
00448 #ifdef TRACE_KSP
00449 block_diag_context->mExtraLAOperations += Timer::GetElapsedTime();
00450 #endif
00451
00452
00453
00454
00455 #ifdef TRACE_KSP
00456 Timer::Reset();
00457 #endif
00458
00459 PetscVecTools::DoInterleavedVecGather(y, block_diag_context->A11_scatter_ctx, block_diag_context->y1_subvector, block_diag_context->A22_scatter_ctx, block_diag_context->y2_subvector);
00460
00461 #ifdef TRACE_KSP
00462 block_diag_context->mGatherTime += Timer::GetElapsedTime();
00463 #endif
00464
00465 return 0;
00466 }