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