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 "PCTwoLevelsBlockDiagonal.hpp"
00037 #include "Exception.hpp"
00038
00039 #include <iostream>
00040
00041 PCTwoLevelsBlockDiagonal::PCTwoLevelsBlockDiagonal(KSP& rKspObject, std::vector<PetscInt>& rBathNodes)
00042 {
00043 PCTwoLevelsBlockDiagonalCreate(rKspObject, rBathNodes);
00044 PCTwoLevelsBlockDiagonalSetUp();
00045 }
00046
00047 PCTwoLevelsBlockDiagonal::~PCTwoLevelsBlockDiagonal()
00048 {
00049 PetscTools::Destroy(mPCContext.A11_matrix_subblock);
00050 PetscTools::Destroy(mPCContext.A22_B1_matrix_subblock);
00051 PetscTools::Destroy(mPCContext.A22_B2_matrix_subblock);
00052
00053 PCDestroy(PETSC_DESTROY_PARAM(mPCContext.PC_amg_A11));
00054 PCDestroy(PETSC_DESTROY_PARAM(mPCContext.PC_amg_A22_B1));
00055 PCDestroy(PETSC_DESTROY_PARAM(mPCContext.PC_amg_A22_B2));
00056
00057 PetscTools::Destroy(mPCContext.x1_subvector);
00058 PetscTools::Destroy(mPCContext.y1_subvector);
00059
00060 PetscTools::Destroy(mPCContext.x21_subvector);
00061 PetscTools::Destroy(mPCContext.y21_subvector);
00062
00063 PetscTools::Destroy(mPCContext.x22_subvector);
00064 PetscTools::Destroy(mPCContext.y22_subvector);
00065
00066 VecScatterDestroy(PETSC_DESTROY_PARAM(mPCContext.A11_scatter_ctx));
00067 VecScatterDestroy(PETSC_DESTROY_PARAM(mPCContext.A22_B1_scatter_ctx));
00068 VecScatterDestroy(PETSC_DESTROY_PARAM(mPCContext.A22_B2_scatter_ctx));
00069 }
00070
00071 void PCTwoLevelsBlockDiagonal::PCTwoLevelsBlockDiagonalCreate(KSP& rKspObject, std::vector<PetscInt>& rBathNodes)
00072 {
00073 KSPGetPC(rKspObject, &mPetscPCObject);
00074
00075 Mat system_matrix, dummy;
00076 MatStructure flag;
00077 KSPGetOperators(rKspObject, &system_matrix, &dummy, &flag);
00078
00079 PetscInt num_rows, num_columns;
00080 MatGetSize(system_matrix, &num_rows, &num_columns);
00081 assert(num_rows==num_columns);
00082
00083 PetscInt num_local_rows, num_local_columns;
00084 MatGetLocalSize(system_matrix, &num_local_rows, &num_local_columns);
00085
00086
00087
00088 if ((num_rows%2 != 0) || (num_local_rows%2 != 0))
00089 {
00090 TERMINATE("Wrong matrix parallel layout detected in PCLDUFactorisation.");
00091 }
00092
00093
00094 unsigned subvector_num_rows = num_rows/2;
00095 unsigned subvector_local_rows = num_local_rows/2;
00096
00097 unsigned subvector_num_rows_tissue = subvector_num_rows - rBathNodes.size();
00098 unsigned subvector_local_rows_tissue = subvector_num_rows_tissue;
00099
00100 unsigned subvector_num_rows_bath = rBathNodes.size();
00101 unsigned subvector_local_rows_bath = subvector_num_rows_bath;
00102
00103 assert(PetscTools::IsSequential());
00104
00105 mPCContext.x1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00106 mPCContext.x21_subvector = PetscTools::CreateVec(subvector_num_rows_tissue, subvector_local_rows_tissue);
00107 mPCContext.x22_subvector = PetscTools::CreateVec(subvector_num_rows_bath, subvector_local_rows_bath);
00108 mPCContext.y1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00109 mPCContext.y21_subvector = PetscTools::CreateVec(subvector_num_rows_tissue, subvector_local_rows_tissue);
00110 mPCContext.y22_subvector = PetscTools::CreateVec(subvector_num_rows_bath, subvector_local_rows_bath);
00111
00112
00113 IS A11_all_rows;
00114 ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 2, &A11_all_rows);
00115
00116 IS A22_all_rows;
00117 PetscInt A22_size;
00118 VecGetSize(mPCContext.x1_subvector, &A22_size);
00120 ISCreateStride(PETSC_COMM_WORLD, A22_size, 1, 2, &A22_all_rows);
00121
00122 IS A22_bath_rows;
00123 PetscInt* phi_e_bath_rows = new PetscInt[rBathNodes.size()];
00124 for (unsigned index=0; index<rBathNodes.size(); index++)
00125 {
00126 phi_e_bath_rows[index] = 2*rBathNodes[index] + 1;
00127 }
00128 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 2) //PETSc 3.2 or later
00129 ISCreateGeneral(PETSC_COMM_WORLD, rBathNodes.size(), phi_e_bath_rows, PETSC_USE_POINTER, &A22_bath_rows);
00130 #else
00131 ISCreateGeneralWithArray(PETSC_COMM_WORLD, rBathNodes.size(), phi_e_bath_rows, &A22_bath_rows);
00132 #endif
00133
00134 IS A22_tissue_rows;
00135 ISDifference(A22_all_rows, A22_bath_rows, &A22_tissue_rows);
00136
00137
00138 {
00139
00140 Vec dummy_vec = PetscTools::CreateVec(num_rows, num_local_rows);
00141
00143 IS& A11_rows=A11_all_rows;
00144 IS& A22_B1_rows=A22_tissue_rows;
00145 IS& A22_B2_rows=A22_bath_rows;
00146
00147 IS all_vector;
00148 ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 1, &all_vector);
00149
00150 IS tissue_vector;
00151 ISCreateStride(PETSC_COMM_WORLD, (num_rows/2)-rBathNodes.size(), 0, 1, &tissue_vector);
00152
00153 IS bath_vector;
00154 ISCreateStride(PETSC_COMM_WORLD, rBathNodes.size(), 0, 1, &bath_vector);
00155
00156 VecScatterCreate(dummy_vec, A11_rows, mPCContext.x1_subvector, all_vector, &mPCContext.A11_scatter_ctx);
00157 VecScatterCreate(dummy_vec, A22_B1_rows, mPCContext.x21_subvector, tissue_vector, &mPCContext.A22_B1_scatter_ctx);
00158 VecScatterCreate(dummy_vec, A22_B2_rows, mPCContext.x22_subvector, bath_vector, &mPCContext.A22_B2_scatter_ctx);
00159
00160 ISDestroy(PETSC_DESTROY_PARAM(all_vector));
00161 ISDestroy(PETSC_DESTROY_PARAM(tissue_vector));
00162 ISDestroy(PETSC_DESTROY_PARAM(bath_vector));
00163
00164 PetscTools::Destroy(dummy_vec);
00165 }
00166
00167
00168 {
00169
00170 PetscInt low, high, global_size;
00171 VecGetOwnershipRange(mPCContext.x1_subvector, &low, &high);
00172 VecGetSize(mPCContext.x1_subvector, &global_size);
00173 assert(global_size == num_rows/2);
00174
00175 IS A11_local_rows;
00176 IS& A11_columns=A11_all_rows;
00177 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &A11_local_rows);
00178
00179 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00180 MatGetSubMatrix(system_matrix, A11_local_rows, A11_columns,
00181 MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00182 #else
00183 MatGetSubMatrix(system_matrix, A11_local_rows, A11_columns, PETSC_DECIDE,
00184 MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00185 #endif
00186
00187 ISDestroy(PETSC_DESTROY_PARAM(A11_local_rows));
00188 }
00189
00190
00191 {
00192
00193
00194
00195
00196
00197
00198 assert(PetscTools::IsSequential());
00199 IS& A22_B1_local_rows = A22_tissue_rows;
00200 IS& A22_B1_columns = A22_tissue_rows;
00201
00202 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00203 MatGetSubMatrix(system_matrix, A22_B1_local_rows, A22_B1_columns,
00204 MAT_INITIAL_MATRIX, &mPCContext.A22_B1_matrix_subblock);
00205 #else
00206 MatGetSubMatrix(system_matrix, A22_B1_local_rows, A22_B1_columns, PETSC_DECIDE,
00207 MAT_INITIAL_MATRIX, &mPCContext.A22_B1_matrix_subblock);
00208 #endif
00209
00210 }
00211
00212
00213 {
00214
00215
00216
00217
00218
00219
00220 assert(PetscTools::IsSequential());
00221 IS& A22_B2_local_rows = A22_bath_rows;
00222 IS& A22_B2_columns = A22_bath_rows;
00223
00224 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00225 MatGetSubMatrix(system_matrix, A22_B2_local_rows, A22_B2_columns,
00226 MAT_INITIAL_MATRIX, &mPCContext.A22_B2_matrix_subblock);
00227 #else
00228 MatGetSubMatrix(system_matrix, A22_B2_local_rows, A22_B2_columns, PETSC_DECIDE,
00229 MAT_INITIAL_MATRIX, &mPCContext.A22_B2_matrix_subblock);
00230 #endif
00231 }
00232
00233 ISDestroy(PETSC_DESTROY_PARAM(A11_all_rows));
00234 ISDestroy(PETSC_DESTROY_PARAM(A22_all_rows));
00235 ISDestroy(PETSC_DESTROY_PARAM(A22_bath_rows));
00236 delete[] phi_e_bath_rows;
00237 ISDestroy(PETSC_DESTROY_PARAM(A22_tissue_rows));
00238
00239
00240 PCSetType(mPetscPCObject, PCSHELL);
00241 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00242 PCShellSetApply(mPetscPCObject, PCTwoLevelsBlockDiagonalApply, (void*) &mPCContext);
00243 #else
00244
00245 PCShellSetContext(mPetscPCObject, &mPCContext);
00246
00247
00248 PCShellSetApply(mPetscPCObject, PCTwoLevelsBlockDiagonalApply);
00249 #endif
00250 }
00251
00252 void PCTwoLevelsBlockDiagonal::PCTwoLevelsBlockDiagonalSetUp()
00253 {
00254
00255 PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00256 PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00257 PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00258
00259
00260 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A11));
00261 PCSetType(mPCContext.PC_amg_A11, PCBJACOBI);
00262 PCSetOperators(mPCContext.PC_amg_A11, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock, DIFFERENT_NONZERO_PATTERN);
00263 PCSetFromOptions(mPCContext.PC_amg_A11);
00264 PCSetUp(mPCContext.PC_amg_A11);
00265
00266
00267 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A22_B1));
00268 PCSetType(mPCContext.PC_amg_A22_B1, PCBJACOBI);
00269 PCSetOperators(mPCContext.PC_amg_A22_B1, mPCContext.A22_B1_matrix_subblock, mPCContext.A22_B1_matrix_subblock, DIFFERENT_NONZERO_PATTERN);
00270 PCSetFromOptions(mPCContext.PC_amg_A22_B1);
00271 PCSetUp(mPCContext.PC_amg_A22_B1);
00272
00273
00274 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A22_B2));
00275 PCSetType(mPCContext.PC_amg_A22_B2, PCHYPRE);
00276
00277 PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00278
00279 PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00280 PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00281 PetscOptionsSetValue("-pc_hypre_boomeramg_coarsen_type", "HMIS");
00282
00283 PCSetOperators(mPCContext.PC_amg_A22_B2, mPCContext.A22_B2_matrix_subblock, mPCContext.A22_B2_matrix_subblock, DIFFERENT_NONZERO_PATTERN);
00284 PCSetFromOptions(mPCContext.PC_amg_A22_B2);
00285 PCSetUp(mPCContext.PC_amg_A22_B2);
00286 }
00287
00288 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
00289 PetscErrorCode PCTwoLevelsBlockDiagonalApply(PC pc_object, Vec x, Vec y)
00290 {
00291 void* pc_context;
00292
00293 PCShellGetContext(pc_object, &pc_context);
00294 #else
00295 PetscErrorCode PCTwoLevelsBlockDiagonalApply(void* pc_context, Vec x, Vec y)
00296 {
00297 #endif
00298
00299
00300 PCTwoLevelsBlockDiagonal::PCTwoLevelsBlockDiagonalContext* block_diag_context = (PCTwoLevelsBlockDiagonal::PCTwoLevelsBlockDiagonalContext*) pc_context;
00301 assert(block_diag_context!=NULL);
00302
00303
00304
00305
00306
00307 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00308 VecScatterBegin(block_diag_context->A11_scatter_ctx, x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD);
00309 VecScatterEnd(block_diag_context->A11_scatter_ctx, x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD);
00310 #else
00311 VecScatterBegin(x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A11_scatter_ctx);
00312 VecScatterEnd(x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A11_scatter_ctx);
00313 #endif
00314
00315
00316 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00317 VecScatterBegin(block_diag_context->A22_B1_scatter_ctx, x, block_diag_context->x21_subvector, INSERT_VALUES, SCATTER_FORWARD);
00318 VecScatterEnd(block_diag_context->A22_B1_scatter_ctx, x, block_diag_context->x21_subvector, INSERT_VALUES, SCATTER_FORWARD);
00319 #else
00320 VecScatterBegin(x, block_diag_context->x21_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A22_B1_scatter_ctx);
00321 VecScatterEnd(x, block_diag_context->x21_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A22_B1_scatter_ctx);
00322 #endif
00323
00324
00325 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00326 VecScatterBegin(block_diag_context->A22_B2_scatter_ctx, x, block_diag_context->x22_subvector, INSERT_VALUES, SCATTER_FORWARD);
00327 VecScatterEnd(block_diag_context->A22_B2_scatter_ctx, x, block_diag_context->x22_subvector, INSERT_VALUES, SCATTER_FORWARD);
00328 #else
00329 VecScatterBegin(x, block_diag_context->x22_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A22_B2_scatter_ctx);
00330 VecScatterEnd(x, block_diag_context->x22_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A22_B2_scatter_ctx);
00331 #endif
00332
00333
00334
00335
00336
00337
00338 PCApply(block_diag_context->PC_amg_A11, block_diag_context->x1_subvector, block_diag_context->y1_subvector);
00339 PCApply(block_diag_context->PC_amg_A22_B1, block_diag_context->x21_subvector, block_diag_context->y21_subvector);
00340 PCApply(block_diag_context->PC_amg_A22_B2, block_diag_context->x22_subvector, block_diag_context->y22_subvector);
00341
00342
00343
00344
00345
00346 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00347 VecScatterBegin(block_diag_context->A11_scatter_ctx, block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00348 VecScatterEnd(block_diag_context->A11_scatter_ctx, block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00349 #else
00350 VecScatterBegin(block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A11_scatter_ctx);
00351 VecScatterEnd(block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A11_scatter_ctx);
00352 #endif
00353
00354
00355 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00356 VecScatterBegin(block_diag_context->A22_B1_scatter_ctx, block_diag_context->y21_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00357 VecScatterEnd(block_diag_context->A22_B1_scatter_ctx, block_diag_context->y21_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00358 #else
00359 VecScatterBegin(block_diag_context->y21_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A22_B1_scatter_ctx);
00360 VecScatterEnd(block_diag_context->y21_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A22_B1_scatter_ctx);
00361 #endif
00362
00363
00364 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00365 VecScatterBegin(block_diag_context->A22_B2_scatter_ctx, block_diag_context->y22_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00366 VecScatterEnd(block_diag_context->A22_B2_scatter_ctx, block_diag_context->y22_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00367 #else
00368 VecScatterBegin(block_diag_context->y22_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A22_B2_scatter_ctx);
00369 VecScatterEnd(block_diag_context->y22_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A22_B2_scatter_ctx);
00370 #endif
00371
00372 return 0;
00373 }