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 "PCBlockDiagonal.hpp"
00040 #include "Exception.hpp"
00041 #include "Warnings.hpp"
00042
00043 PCBlockDiagonal::PCBlockDiagonal(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 PCBlockDiagonalCreate(rKspObject);
00053 PCBlockDiagonalSetUp();
00054 }
00055
00056 PCBlockDiagonal::~PCBlockDiagonal()
00057 {
00058 #ifdef TRACE_KSP
00059 if (PetscTools::AmMaster())
00060 {
00061 std::cout << " -- Block diagonal 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 mGatherTime: " << mPCContext.mGatherTime << std::endl;
00066 }
00067 #endif
00068
00069 PetscTools::Destroy(mPCContext.A11_matrix_subblock);
00070 PetscTools::Destroy(mPCContext.A22_matrix_subblock);
00071
00072 PCDestroy(PETSC_DESTROY_PARAM(mPCContext.PC_amg_A11));
00073 PCDestroy(PETSC_DESTROY_PARAM(mPCContext.PC_amg_A22));
00074
00075 PetscTools::Destroy(mPCContext.x1_subvector);
00076 PetscTools::Destroy(mPCContext.y1_subvector);
00077
00078 PetscTools::Destroy(mPCContext.x2_subvector);
00079 PetscTools::Destroy(mPCContext.y2_subvector);
00080
00081 VecScatterDestroy(PETSC_DESTROY_PARAM(mPCContext.A11_scatter_ctx));
00082 VecScatterDestroy(PETSC_DESTROY_PARAM(mPCContext.A22_scatter_ctx));
00083 }
00084
00085 void PCBlockDiagonal::PCBlockDiagonalCreate(KSP& rKspObject)
00086 {
00087 KSPGetPC(rKspObject, &mPetscPCObject);
00088
00089 Mat system_matrix, dummy;
00090 MatStructure flag;
00091 KSPGetOperators(rKspObject, &system_matrix, &dummy, &flag);
00092
00093 PetscInt num_rows, num_columns;
00094 MatGetSize(system_matrix, &num_rows, &num_columns);
00095 assert(num_rows==num_columns);
00096
00097 PetscInt num_local_rows, num_local_columns;
00098 MatGetLocalSize(system_matrix, &num_local_rows, &num_local_columns);
00099
00100
00101
00102 if ((num_rows%2 != 0) || (num_local_rows%2 != 0))
00103 {
00104 TERMINATE("Wrong matrix parallel layout detected in PCLDUFactorisation.");
00105 }
00106
00107
00108 unsigned subvector_num_rows = num_rows/2;
00109 unsigned subvector_local_rows = num_local_rows/2;
00110 mPCContext.x1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00111 mPCContext.x2_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00112 mPCContext.y1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00113 mPCContext.y2_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00114
00115
00116 {
00117
00118 Vec dummy_vec = PetscTools::CreateVec(num_rows, num_local_rows);
00119
00120 PetscVecTools::SetupInterleavedVectorScatterGather(dummy_vec, mPCContext.A11_scatter_ctx, mPCContext.A22_scatter_ctx);
00121
00122 PetscTools::Destroy(dummy_vec);
00123 }
00124
00125
00126 {
00127
00128 PetscInt low, high, global_size;
00129 VecGetOwnershipRange(mPCContext.x1_subvector, &low, &high);
00130 VecGetSize(mPCContext.x1_subvector, &global_size);
00131 assert(global_size == num_rows/2);
00132
00133 IS A11_local_rows;
00134 IS A11_columns;
00135 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &A11_local_rows);
00136 ISCreateStride(PETSC_COMM_WORLD, global_size, 0, 2, &A11_columns);
00137
00138 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00139 MatGetSubMatrix(system_matrix, A11_local_rows, A11_local_rows,
00140 MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00141 #else
00142 MatGetSubMatrix(system_matrix, A11_local_rows, A11_columns, PETSC_DECIDE,
00143 MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00144 #endif
00145
00146
00147 ISDestroy(PETSC_DESTROY_PARAM(A11_local_rows));
00148 ISDestroy(PETSC_DESTROY_PARAM(A11_columns));
00149 }
00150
00151
00152 {
00153
00154 PetscInt low, high, global_size;
00155 VecGetOwnershipRange(mPCContext.x2_subvector, &low, &high);
00156 VecGetSize(mPCContext.x2_subvector, &global_size);
00157 assert(global_size == num_rows/2);
00158
00159 IS A22_local_rows;
00160 IS A22_columns;
00161 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low+1, 2, &A22_local_rows);
00162 ISCreateStride(PETSC_COMM_WORLD, global_size, 1, 2, &A22_columns);
00163
00164 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00165 MatGetSubMatrix(system_matrix, A22_local_rows, A22_local_rows,
00166 MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00167 #else
00168 MatGetSubMatrix(system_matrix, A22_local_rows, A22_columns, PETSC_DECIDE,
00169 MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00170 #endif
00171
00172 ISDestroy(PETSC_DESTROY_PARAM(A22_local_rows));
00173 ISDestroy(PETSC_DESTROY_PARAM(A22_columns));
00174 }
00175
00176
00177 PCSetType(mPetscPCObject, PCSHELL);
00178 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00179 PCShellSetApply(mPetscPCObject, PCBlockDiagonalApply, (void*) &mPCContext);
00180 #else
00181
00182 PCShellSetContext(mPetscPCObject, &mPCContext);
00183
00184
00185 PCShellSetApply(mPetscPCObject, PCBlockDiagonalApply);
00186 #endif
00187
00188 }
00189
00190 void PCBlockDiagonal::PCBlockDiagonalSetUp()
00191 {
00192
00193
00194
00195
00196
00197
00198 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A11));
00199
00200
00201
00202
00203
00204 PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL);
00205 PetscErrorCode pc_set_error = PCSetType(mPCContext.PC_amg_A11, PCHYPRE);
00206 if (pc_set_error != 0)
00207 {
00208 WARNING("PETSc hypre preconditioning library is not installed");
00209 }
00210
00211 PetscPopErrorHandler();
00212
00213
00214
00215 PetscOptionsSetValue("-pc_hypre_type", "euclid");
00216 PetscOptionsSetValue("-pc_hypre_euclid_levels", "0");
00217
00218
00219
00220
00221
00222
00223
00224
00225 PCSetOperators(mPCContext.PC_amg_A11, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock, SAME_PRECONDITIONER);
00226 PCSetFromOptions(mPCContext.PC_amg_A11);
00227 PCSetUp(mPCContext.PC_amg_A11);
00228
00229
00230 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A22));
00231
00232
00233 PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL);
00234 PCSetType(mPCContext.PC_amg_A22, PCHYPRE);
00235
00236 PetscPopErrorHandler();
00237
00238
00239 PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00240 PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00241 PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00242 PetscOptionsSetValue("-pc_hypre_boomeramg_coarsen_type", "HMIS");
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272
00273
00274
00275
00276
00277
00278
00279 PCSetOperators(mPCContext.PC_amg_A22, mPCContext.A22_matrix_subblock, mPCContext.A22_matrix_subblock, SAME_PRECONDITIONER);
00280 PCSetFromOptions(mPCContext.PC_amg_A22);
00281 PCSetUp(mPCContext.PC_amg_A22);
00282 }
00283
00284 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00285 PetscErrorCode PCBlockDiagonalApply(PC pc_object, Vec x, Vec y)
00286 {
00287 void* pc_context;
00288
00289 PCShellGetContext(pc_object, &pc_context);
00290 #else
00291 PetscErrorCode PCBlockDiagonalApply(void* pc_context, Vec x, Vec y)
00292 {
00293 #endif
00294
00295
00296 PCBlockDiagonal::PCBlockDiagonalContext* block_diag_context = (PCBlockDiagonal::PCBlockDiagonalContext*) pc_context;
00297 assert(block_diag_context!=NULL);
00298
00299
00300
00301
00302 #ifdef TRACE_KSP
00303 double init_time = MPI_Wtime();
00304 #endif
00305
00306 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);
00307
00308 #ifdef TRACE_KSP
00309 block_diag_context->mScatterTime += MPI_Wtime() - init_time;
00310 #endif
00311
00312
00313
00314
00315
00316 #ifdef TRACE_KSP
00317 init_time = MPI_Wtime();
00318 #endif
00319 PCApply(block_diag_context->PC_amg_A11, block_diag_context->x1_subvector, block_diag_context->y1_subvector);
00320 #ifdef TRACE_KSP
00321 block_diag_context->mA1PreconditionerTime += MPI_Wtime() - init_time;
00322 #endif
00323
00324 #ifdef TRACE_KSP
00325 init_time = MPI_Wtime();
00326 #endif
00327 PCApply(block_diag_context->PC_amg_A22, block_diag_context->x2_subvector, block_diag_context->y2_subvector);
00328 #ifdef TRACE_KSP
00329 block_diag_context->mA2PreconditionerTime += MPI_Wtime() - init_time;
00330 #endif
00331
00332
00333
00334
00335
00336 #ifdef TRACE_KSP
00337 init_time = MPI_Wtime();
00338 #endif
00339
00340 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);
00341
00342 #ifdef TRACE_KSP
00343 block_diag_context->mGatherTime += MPI_Wtime() - init_time;
00344 #endif
00345 return 0;
00346 }