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 #include "PCLDUFactorisation.hpp"
00030 #include "Exception.hpp"
00031
00032 PCLDUFactorisation::PCLDUFactorisation(KSP& rKspObject)
00033 {
00034 PCLDUFactorisationCreate(rKspObject);
00035 PCLDUFactorisationSetUp();
00036 }
00037
00038 PCLDUFactorisation::~PCLDUFactorisation()
00039 {
00040 MatDestroy(mPCContext.A11_matrix_subblock);
00041 MatDestroy(mPCContext.A22_matrix_subblock);
00042 MatDestroy(mPCContext.B_matrix_subblock);
00043
00044 PCDestroy(mPCContext.PC_amg_A11);
00045 PCDestroy(mPCContext.PC_amg_A22);
00046
00047 VecDestroy(mPCContext.x1_subvector);
00048 VecDestroy(mPCContext.y1_subvector);
00049
00050 VecDestroy(mPCContext.x2_subvector);
00051 VecDestroy(mPCContext.y2_subvector);
00052
00053 VecDestroy(mPCContext.z);
00054 VecDestroy(mPCContext.temp);
00055 }
00056
00057 void PCLDUFactorisation::PCLDUFactorisationCreate(KSP& rKspObject)
00058 {
00059 KSPGetPC(rKspObject, &mPetscPCObject);
00060
00061 Mat system_matrix, dummy;
00062 MatStructure flag;
00063 KSPGetOperators(rKspObject, &system_matrix, &dummy, &flag);
00064
00065 PetscInt num_rows, num_columns;
00066 MatGetSize(system_matrix, &num_rows, &num_columns);
00067 assert(num_rows==num_columns);
00068
00069 PCSetType(mPetscPCObject, PCSHELL);
00070 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00071
00072 PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply, (void*) &mPCContext);
00073 #else
00074
00075 PCShellSetContext(mPetscPCObject, &mPCContext);
00076
00077 PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply);
00078 #endif
00079
00080
00081 IS A11_rows, A11_columns;
00082 ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 2, &A11_rows);
00083 ISCreateStride(PETSC_COMM_WORLD, num_columns/2, 0, 2, &A11_columns);
00084
00085 MatGetSubMatrix(system_matrix, A11_rows, A11_columns, PETSC_DECIDE, MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00086
00087 ISDestroy(A11_rows);
00088 ISDestroy(A11_columns);
00089
00090
00091 IS A22_rows, A22_columns;
00092 ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 1, 2, &A22_rows);
00093 ISCreateStride(PETSC_COMM_WORLD, num_columns/2, 1, 2, &A22_columns);
00094
00095 MatGetSubMatrix(system_matrix, A22_rows, A22_columns, PETSC_DECIDE, MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00096
00097 ISDestroy(A22_rows);
00098 ISDestroy(A22_columns);
00099
00100
00101 IS B_rows, B_columns;
00102 ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 2, &B_rows);
00103 ISCreateStride(PETSC_COMM_WORLD, num_columns/2, 1, 2, &B_columns);
00104
00105 MatGetSubMatrix(system_matrix, B_rows, B_columns, PETSC_DECIDE, MAT_INITIAL_MATRIX, &mPCContext.B_matrix_subblock);
00106
00107 ISDestroy(B_rows);
00108 ISDestroy(B_columns);
00109
00110
00111 mPCContext.x1_subvector = PetscTools::CreateVec(num_rows/2);
00112 mPCContext.x2_subvector = PetscTools::CreateVec(num_rows/2);
00113 mPCContext.y1_subvector = PetscTools::CreateVec(num_rows/2);
00114 mPCContext.y2_subvector = PetscTools::CreateVec(num_rows/2);
00115 mPCContext.z = PetscTools::CreateVec(num_rows/2);
00116 mPCContext.temp = PetscTools::CreateVec(num_rows/2);
00117 }
00118
00119 void PCLDUFactorisation::PCLDUFactorisationSetUp()
00120 {
00121
00122 PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00123 PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00124 PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00125
00126
00127
00128
00129 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A11));
00130 PCSetOperators(mPCContext.PC_amg_A11, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock, DIFFERENT_NONZERO_PATTERN);
00131
00132
00133
00135 PCSetType(mPCContext.PC_amg_A11, PCHYPRE);
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00153
00154 PCSetFromOptions(mPCContext.PC_amg_A11);
00155 PCSetUp(mPCContext.PC_amg_A11);
00156
00157
00158
00159
00160
00161 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A22));
00162 PCSetOperators(mPCContext.PC_amg_A22, mPCContext.A22_matrix_subblock, mPCContext.A22_matrix_subblock, DIFFERENT_NONZERO_PATTERN);
00163
00164
00165
00167 PCSetType(mPCContext.PC_amg_A22, PCHYPRE);
00169
00170
00171
00172
00173
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183
00185
00186 PCSetFromOptions(mPCContext.PC_amg_A22);
00187 PCSetUp(mPCContext.PC_amg_A22);
00188 }
00189
00190 PetscErrorCode PCLDUFactorisationApply(void* pc_context, Vec x, Vec y)
00191 {
00193
00194
00195 PCLDUFactorisation::PCLDUFactorisationContext* block_diag_context = (PCLDUFactorisation::PCLDUFactorisationContext*) pc_context;
00196 assert(block_diag_context!=NULL);
00197
00198
00199
00200
00201 PetscInt num_rows;
00202 VecGetSize(x, &num_rows);
00203
00204 IS A11_rows;
00205 ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 2, &A11_rows);
00206
00207 VecScatter A11_scatter_ctx;
00208 VecScatterCreate(x, A11_rows, block_diag_context->x1_subvector, PETSC_NULL, &A11_scatter_ctx);
00209
00210
00211 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00212 VecScatterBegin(A11_scatter_ctx, x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD);
00213 VecScatterEnd(A11_scatter_ctx, x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD);
00214 #else
00215 VecScatterBegin(x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD, A11_scatter_ctx);
00216 VecScatterEnd(x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD, A11_scatter_ctx);
00217 #endif
00218
00219 IS A22_rows;
00220 ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 1, 2, &A22_rows);
00221
00222 VecScatter A22_scatter_ctx;
00223 VecScatterCreate(x, A22_rows, block_diag_context->x2_subvector, PETSC_NULL, &A22_scatter_ctx);
00224
00225
00226 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00227 VecScatterBegin(A22_scatter_ctx, x, block_diag_context->x2_subvector, INSERT_VALUES, SCATTER_FORWARD);
00228 VecScatterEnd(A22_scatter_ctx, x, block_diag_context->x2_subvector, INSERT_VALUES, SCATTER_FORWARD);
00229 #else
00230 VecScatterBegin(x, block_diag_context->x2_subvector, INSERT_VALUES, SCATTER_FORWARD, A22_scatter_ctx);
00231 VecScatterEnd(x, block_diag_context->x2_subvector, INSERT_VALUES, SCATTER_FORWARD, A22_scatter_ctx);
00232 #endif
00233
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243 PCApply(block_diag_context->PC_amg_A11, block_diag_context->x1_subvector, block_diag_context->z);
00244
00245
00246 MatMult(block_diag_context->B_matrix_subblock,block_diag_context->z,block_diag_context->temp);
00247 double minus_one = -1.0;
00248 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00249 VecAYPX(&minus_one, block_diag_context->x2_subvector, block_diag_context->temp);
00250 #else
00251 VecAYPX(block_diag_context->temp, minus_one, block_diag_context->x2_subvector);
00252 #endif
00253 PCApply(block_diag_context->PC_amg_A22, block_diag_context->temp, block_diag_context->y2_subvector);
00254
00255
00256 MatMult(block_diag_context->B_matrix_subblock,block_diag_context->y2_subvector,block_diag_context->temp);
00257 PCApply(block_diag_context->PC_amg_A11, block_diag_context->temp, block_diag_context->y1_subvector);
00258 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00259 VecAYPX(&minus_one, block_diag_context->z, block_diag_context->y1_subvector);
00260 #else
00261 VecAYPX(block_diag_context->y1_subvector, minus_one, block_diag_context->z);
00262 #endif
00263
00264
00265
00266
00267
00268
00269 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00270 VecScatterBegin(A11_scatter_ctx, block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00271 VecScatterEnd(A11_scatter_ctx, block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00272 #else
00273 VecScatterBegin(block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE, A11_scatter_ctx);
00274 VecScatterEnd(block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE, A11_scatter_ctx);
00275 #endif
00276
00277
00278 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00279 VecScatterBegin(A22_scatter_ctx, block_diag_context->y2_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00280 VecScatterEnd(A22_scatter_ctx, block_diag_context->y2_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00281 #else
00282 VecScatterBegin(block_diag_context->y2_subvector, y, INSERT_VALUES, SCATTER_REVERSE, A22_scatter_ctx);
00283 VecScatterEnd(block_diag_context->y2_subvector, y, INSERT_VALUES, SCATTER_REVERSE, A22_scatter_ctx);
00284 #endif
00285
00286
00287
00288
00289 ISDestroy(A11_rows);
00290 ISDestroy(A22_rows);
00291
00292 VecScatterDestroy(A11_scatter_ctx);
00293 VecScatterDestroy(A22_scatter_ctx);
00294
00295 return 0;
00296 }