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 MatStructure flag;
00094 KSPGetOperators(rKspObject, &system_matrix, &dummy, &flag);
00095
00096 PetscInt num_rows, num_columns;
00097 MatGetSize(system_matrix, &num_rows, &num_columns);
00098
00099 PetscInt num_local_rows, num_local_columns;
00100 MatGetLocalSize(system_matrix, &num_local_rows, &num_local_columns);
00101
00102
00103
00104 if ((num_rows%2 != 0) || (num_local_rows%2 != 0))
00105 {
00106 TERMINATE("Wrong matrix parallel layout detected in PCLDUFactorisation.");
00107 }
00108
00109
00110 unsigned subvector_num_rows = num_rows/2;
00111 unsigned subvector_local_rows = num_local_rows/2;
00112 mPCContext.x1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00113 mPCContext.x2_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00114 mPCContext.y1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00115 mPCContext.y2_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00116 mPCContext.z = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00117 mPCContext.temp = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00118
00119
00120 {
00121
00122 Vec dummy_vec = PetscTools::CreateVec(num_rows, num_local_rows);
00123
00124 PetscVecTools::SetupInterleavedVectorScatterGather(dummy_vec, mPCContext.A11_scatter_ctx, mPCContext.A22_scatter_ctx);
00125
00126 PetscTools::Destroy(dummy_vec);
00127 }
00128
00129
00130 {
00131
00132 PetscInt low, high, global_size;
00133 VecGetOwnershipRange(mPCContext.x1_subvector, &low, &high);
00134 VecGetSize(mPCContext.x1_subvector, &global_size);
00135 assert(global_size == num_rows/2);
00136
00137 IS A11_local_rows;
00138 IS A11_columns;
00139 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &A11_local_rows);
00140 ISCreateStride(PETSC_COMM_WORLD, global_size, 0, 2, &A11_columns);
00141
00142 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00143 MatGetSubMatrix(system_matrix, A11_local_rows, A11_local_rows,
00144 MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00145 #else
00146 MatGetSubMatrix(system_matrix, A11_local_rows, A11_columns, PETSC_DECIDE,
00147 MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00148 #endif
00149 ISDestroy(PETSC_DESTROY_PARAM(A11_local_rows));
00150 ISDestroy(PETSC_DESTROY_PARAM(A11_columns));
00151 }
00152
00153
00154 {
00155
00156 PetscInt low, high, global_size;
00157 VecGetOwnershipRange(mPCContext.x2_subvector, &low, &high);
00158 VecGetSize(mPCContext.x2_subvector, &global_size);
00159 assert(global_size == num_rows/2);
00160
00161 IS A22_local_rows;
00162 IS A22_columns;
00163 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low+1, 2, &A22_local_rows);
00164 ISCreateStride(PETSC_COMM_WORLD, global_size, 1, 2, &A22_columns);
00165
00166 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00167 MatGetSubMatrix(system_matrix, A22_local_rows, A22_local_rows,
00168 MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00169 #else
00170 MatGetSubMatrix(system_matrix, A22_local_rows, A22_columns, PETSC_DECIDE,
00171 MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00172 #endif
00173
00174 ISDestroy(PETSC_DESTROY_PARAM(A22_local_rows));
00175 ISDestroy(PETSC_DESTROY_PARAM(A22_columns));
00176 }
00177
00178
00179 {
00180
00181 PetscInt low, high, global_size;
00182 VecGetOwnershipRange(mPCContext.x1_subvector, &low, &high);
00183 VecGetSize(mPCContext.x1_subvector, &global_size);
00184 assert(global_size == num_rows/2);
00185
00186 IS B_local_rows;
00187 IS B_columns;
00188 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &B_local_rows);
00189
00190 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00191 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low+1, 2, &B_columns);
00192 MatGetSubMatrix(system_matrix, B_local_rows, B_columns,
00193 MAT_INITIAL_MATRIX, &mPCContext.B_matrix_subblock);
00194 #else
00195 ISCreateStride(PETSC_COMM_WORLD, global_size, 1, 2, &B_columns);
00196 MatGetSubMatrix(system_matrix, B_local_rows, B_columns, PETSC_DECIDE,
00197 MAT_INITIAL_MATRIX, &mPCContext.B_matrix_subblock);
00198 #endif
00199
00200 ISDestroy(PETSC_DESTROY_PARAM(B_local_rows));
00201 ISDestroy(PETSC_DESTROY_PARAM(B_columns));
00202 }
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219
00220
00221
00222
00223
00224 PCSetType(mPetscPCObject, PCSHELL);
00225 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00226
00227 PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply, (void*) &mPCContext);
00228 #else
00229
00230 PCShellSetContext(mPetscPCObject, &mPCContext);
00231
00232 PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply);
00233 #endif
00234 }
00235
00236 void PCLDUFactorisation::PCLDUFactorisationSetUp()
00237 {
00238
00239
00240
00241
00242
00243
00244
00245
00246 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A11));
00247 PCSetOperators(mPCContext.PC_amg_A11, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock, SAME_PRECONDITIONER);
00248
00249
00250
00252
00253
00254
00255
00256 PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL);
00257 PetscErrorCode pc_set_error = PCSetType(mPCContext.PC_amg_A11, PCHYPRE);
00258 if (pc_set_error != 0)
00259 {
00260 WARNING("PETSc hypre preconditioning library is not installed");
00261 }
00262
00263 PetscPopErrorHandler();
00264
00265
00266 PetscOptionsSetValue("-pc_hypre_type", "euclid");
00267 PetscOptionsSetValue("-pc_hypre_euclid_levels", "0");
00268
00269
00270
00271
00272
00273
00274
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00292
00293 PCSetFromOptions(mPCContext.PC_amg_A11);
00294 PCSetUp(mPCContext.PC_amg_A11);
00295
00296
00297
00298
00299 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A22));
00300 PCSetOperators(mPCContext.PC_amg_A22, mPCContext.A22_matrix_subblock, mPCContext.A22_matrix_subblock, SAME_PRECONDITIONER);
00301
00302
00303
00305
00306
00307 PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL);
00308 PCSetType(mPCContext.PC_amg_A22, PCHYPRE);
00309
00310 PetscPopErrorHandler();
00311
00312 PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00313 PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00314 PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00315 PetscOptionsSetValue("-pc_hypre_boomeramg_coarsen_type", "HMIS");
00316
00317
00319
00320
00321
00322
00323
00324
00325
00326
00327
00328
00329
00330
00331
00332
00333
00335
00336 PCSetFromOptions(mPCContext.PC_amg_A22);
00337 PCSetUp(mPCContext.PC_amg_A22);
00338 }
00339 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00340 PetscErrorCode PCLDUFactorisationApply(PC pc_object, Vec x, Vec y)
00341 {
00342 void* pc_context;
00343
00344 PCShellGetContext(pc_object, &pc_context);
00345 #else
00346 PetscErrorCode PCLDUFactorisationApply(void* pc_context, Vec x, Vec y)
00347 {
00348 #endif
00350
00351
00352 PCLDUFactorisation::PCLDUFactorisationContext* block_diag_context = (PCLDUFactorisation::PCLDUFactorisationContext*) pc_context;
00353 assert(block_diag_context!=NULL);
00354
00355
00356
00357
00358 #ifdef TRACE_KSP
00359 double init_time = MPI_Wtime();
00360 #endif
00361
00362 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);
00363
00364 #ifdef TRACE_KSP
00365 block_diag_context->mScatterTime += MPI_Wtime() - init_time;
00366 #endif
00367
00368
00369
00370
00371
00372
00373
00374
00375 #ifdef TRACE_KSP
00376 init_time = MPI_Wtime();
00377 #endif
00378
00379 PCApply(block_diag_context->PC_amg_A11, block_diag_context->x1_subvector, block_diag_context->z);
00380 #ifdef TRACE_KSP
00381 block_diag_context->mA1PreconditionerTime += MPI_Wtime() - init_time;
00382 #endif
00383
00384
00385 #ifdef TRACE_KSP
00386 init_time = MPI_Wtime();
00387 #endif
00388 MatMult(block_diag_context->B_matrix_subblock,block_diag_context->z,block_diag_context->temp);
00389 double minus_one = -1.0;
00390 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00391 VecAYPX(&minus_one, block_diag_context->x2_subvector, block_diag_context->temp);
00392 #else
00393 VecAYPX(block_diag_context->temp, minus_one, block_diag_context->x2_subvector);
00394 #endif
00395 #ifdef TRACE_KSP
00396 block_diag_context->mExtraLAOperations += MPI_Wtime() - init_time;
00397 #endif
00398
00399
00400 #ifdef TRACE_KSP
00401 init_time = MPI_Wtime();
00402 #endif
00403 PCApply(block_diag_context->PC_amg_A22, block_diag_context->temp, block_diag_context->y2_subvector);
00404 #ifdef TRACE_KSP
00405 block_diag_context->mA2PreconditionerTime += MPI_Wtime() - init_time;
00406 #endif
00407
00408
00409 #ifdef TRACE_KSP
00410 init_time = MPI_Wtime();
00411 #endif
00412 MatMult(block_diag_context->B_matrix_subblock,block_diag_context->y2_subvector,block_diag_context->temp);
00413 #ifdef TRACE_KSP
00414 block_diag_context->mExtraLAOperations += MPI_Wtime() - init_time;
00415 #endif
00416 #ifdef TRACE_KSP
00417 init_time = MPI_Wtime();
00418 #endif
00419 PCApply(block_diag_context->PC_amg_A11, block_diag_context->temp, block_diag_context->y1_subvector);
00420 #ifdef TRACE_KSP
00421 block_diag_context->mA1PreconditionerTime += MPI_Wtime() - init_time;
00422 #endif
00423
00424 #ifdef TRACE_KSP
00425 init_time = MPI_Wtime();
00426 #endif
00427 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00428 VecAYPX(&minus_one, block_diag_context->z, block_diag_context->y1_subvector);
00429 #else
00430 VecAYPX(block_diag_context->y1_subvector, minus_one, block_diag_context->z);
00431 #endif
00432 #ifdef TRACE_KSP
00433 block_diag_context->mExtraLAOperations += MPI_Wtime() - init_time;
00434 #endif
00435
00436
00437
00438
00439 #ifdef TRACE_KSP
00440 init_time = MPI_Wtime();
00441 #endif
00442
00443 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);
00444
00445 #ifdef TRACE_KSP
00446 block_diag_context->mGatherTime += MPI_Wtime() - init_time;
00447 #endif
00448
00449 return 0;
00450 }