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 <iostream>
00030
00031 #include "PetscVecTools.hpp"
00032 #include "PCLDUFactorisation.hpp"
00033 #include "Exception.hpp"
00034
00035 PCLDUFactorisation::PCLDUFactorisation(KSP& rKspObject)
00036 {
00037 #ifdef TRACE_KSP
00038 mPCContext.mScatterTime = 0.0;
00039 mPCContext.mA1PreconditionerTime = 0.0;;
00040 mPCContext.mA2PreconditionerTime = 0.0;;
00041 mPCContext.mGatherTime = 0.0;;
00042 #endif
00043
00044 PCLDUFactorisationCreate(rKspObject);
00045 PCLDUFactorisationSetUp();
00046 }
00047
00048 PCLDUFactorisation::~PCLDUFactorisation()
00049 {
00050 #ifdef TRACE_KSP
00051 if (PetscTools::AmMaster())
00052 {
00053 std::cout << " -- LDU factorisation preconditioner profile information: " << std::endl;
00054 std::cout << "\t mScatterTime: " << mPCContext.mScatterTime << std::endl;
00055 std::cout << "\t mA1PreconditionerTime: " << mPCContext.mA1PreconditionerTime << std::endl;
00056 std::cout << "\t mA2PreconditionerTime: " << mPCContext.mA2PreconditionerTime << std::endl;
00057 std::cout << "\t mExtraLAOperations: " << mPCContext.mExtraLAOperations << std::endl;
00058 std::cout << "\t mGatherTime: " << mPCContext.mGatherTime << std::endl;
00059 }
00060 #endif
00061
00062 MatDestroy(mPCContext.A11_matrix_subblock);
00063 MatDestroy(mPCContext.A22_matrix_subblock);
00064 MatDestroy(mPCContext.B_matrix_subblock);
00065
00066 PCDestroy(mPCContext.PC_amg_A11);
00067 PCDestroy(mPCContext.PC_amg_A22);
00068
00069 VecDestroy(mPCContext.x1_subvector);
00070 VecDestroy(mPCContext.y1_subvector);
00071 VecDestroy(mPCContext.x2_subvector);
00072 VecDestroy(mPCContext.y2_subvector);
00073 VecDestroy(mPCContext.z);
00074 VecDestroy(mPCContext.temp);
00075
00076 VecScatterDestroy(mPCContext.A11_scatter_ctx);
00077 VecScatterDestroy(mPCContext.A22_scatter_ctx);
00078 }
00079
00080 void PCLDUFactorisation::PCLDUFactorisationCreate(KSP& rKspObject)
00081 {
00082 KSPGetPC(rKspObject, &mPetscPCObject);
00083
00084 Mat system_matrix, dummy;
00085 MatStructure flag;
00086 KSPGetOperators(rKspObject, &system_matrix, &dummy, &flag);
00087
00088 PetscInt num_rows, num_columns;
00089 MatGetSize(system_matrix, &num_rows, &num_columns);
00090
00091 PetscInt num_local_rows, num_local_columns;
00092 MatGetLocalSize(system_matrix, &num_local_rows, &num_local_columns);
00093
00094
00095
00096 if ((num_rows%2 != 0) || (num_local_rows%2 != 0))
00097 {
00098 TERMINATE("Wrong matrix parallel layout detected in PCLDUFactorisation.");
00099 }
00100
00101
00102 unsigned subvector_num_rows = num_rows/2;
00103 unsigned subvector_local_rows = num_local_rows/2;
00104 mPCContext.x1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00105 mPCContext.x2_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00106 mPCContext.y1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00107 mPCContext.y2_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00108 mPCContext.z = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00109 mPCContext.temp = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00110
00111
00112 {
00113
00114 Vec dummy_vec = PetscTools::CreateVec(num_rows, num_local_rows);
00115
00116 PetscVecTools::SetupInterleavedVectorScatterGather(dummy_vec, mPCContext.A11_scatter_ctx, mPCContext.A22_scatter_ctx);
00117
00118 VecDestroy(dummy_vec);
00119 }
00120
00121
00122
00123 {
00124
00125 PetscInt low, high, global_size;
00126 VecGetOwnershipRange(mPCContext.x1_subvector, &low, &high);
00127 VecGetSize(mPCContext.x1_subvector, &global_size);
00128 assert(global_size == num_rows/2);
00129
00130 IS A11_local_rows;
00131 IS A11_columns;
00132 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &A11_local_rows);
00133 ISCreateStride(PETSC_COMM_WORLD, global_size, 0, 2, &A11_columns);
00134
00135 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00136 MatGetSubMatrix(system_matrix, A11_local_rows, A11_columns,
00137 MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00138 #else
00139 MatGetSubMatrix(system_matrix, A11_local_rows, A11_columns, PETSC_DECIDE,
00140 MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00141 #endif
00142 ISDestroy(A11_local_rows);
00143 ISDestroy(A11_columns);
00144 }
00145
00146
00147 {
00148
00149 PetscInt low, high, global_size;
00150 VecGetOwnershipRange(mPCContext.x2_subvector, &low, &high);
00151 VecGetSize(mPCContext.x2_subvector, &global_size);
00152 assert(global_size == num_rows/2);
00153
00154 IS A22_local_rows;
00155 IS A22_columns;
00156 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low+1, 2, &A22_local_rows);
00157 ISCreateStride(PETSC_COMM_WORLD, global_size, 1, 2, &A22_columns);
00158
00159 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00160 MatGetSubMatrix(system_matrix, A22_local_rows, A22_columns,
00161 MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00162 #else
00163 MatGetSubMatrix(system_matrix, A22_local_rows, A22_columns, PETSC_DECIDE,
00164 MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00165 #endif
00166
00167 ISDestroy(A22_local_rows);
00168 ISDestroy(A22_columns);
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 B_local_rows;
00180 IS B_columns;
00181 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &B_local_rows);
00182 ISCreateStride(PETSC_COMM_WORLD, global_size, 1, 2, &B_columns);
00183
00184 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00185 MatGetSubMatrix(system_matrix, B_local_rows, B_columns,
00186 MAT_INITIAL_MATRIX, &mPCContext.B_matrix_subblock);
00187 #else
00188 MatGetSubMatrix(system_matrix, B_local_rows, B_columns, PETSC_DECIDE,
00189 MAT_INITIAL_MATRIX, &mPCContext.B_matrix_subblock);
00190 #endif
00191
00192 ISDestroy(B_local_rows);
00193 ISDestroy(B_columns);
00194 }
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218 PCSetType(mPetscPCObject, PCSHELL);
00219 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00220
00221 PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply, (void*) &mPCContext);
00222 #else
00223
00224 PCShellSetContext(mPetscPCObject, &mPCContext);
00225
00226 PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply);
00227 #endif
00228
00229 }
00230
00231 void PCLDUFactorisation::PCLDUFactorisationSetUp()
00232 {
00233
00234
00235
00236
00237
00238
00239
00240
00241 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A11));
00242 PCSetOperators(mPCContext.PC_amg_A11, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock, SAME_PRECONDITIONER);
00243
00244
00245
00247
00248
00249 PCSetType(mPCContext.PC_amg_A11, PCHYPRE);
00250
00251 PetscOptionsSetValue("-pc_hypre_type", "euclid");
00252 PetscOptionsSetValue("-pc_hypre_euclid_levels", "0");
00253
00254
00255
00256
00257
00258
00259
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272
00273
00274
00275
00277
00278 PCSetFromOptions(mPCContext.PC_amg_A11);
00279 PCSetUp(mPCContext.PC_amg_A11);
00280
00281
00282
00283
00284
00285 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A22));
00286 PCSetOperators(mPCContext.PC_amg_A22, mPCContext.A22_matrix_subblock, mPCContext.A22_matrix_subblock, SAME_PRECONDITIONER);
00287
00288
00289
00291 PCSetType(mPCContext.PC_amg_A22, PCHYPRE);
00292 PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00293 PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00294 PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00295 PetscOptionsSetValue("-pc_hypre_boomeramg_coarsen_type", "HMIS");
00296
00297
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00315
00316 PCSetFromOptions(mPCContext.PC_amg_A22);
00317 PCSetUp(mPCContext.PC_amg_A22);
00318 }
00319 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00320 PetscErrorCode PCLDUFactorisationApply(PC pc_object, Vec x, Vec y)
00321 {
00322 void* pc_context;
00323
00324 PCShellGetContext(pc_object, &pc_context);
00325 #else
00326 PetscErrorCode PCLDUFactorisationApply(void* pc_context, Vec x, Vec y)
00327 {
00328 #endif
00330
00331
00332 PCLDUFactorisation::PCLDUFactorisationContext* block_diag_context = (PCLDUFactorisation::PCLDUFactorisationContext*) pc_context;
00333 assert(block_diag_context!=NULL);
00334
00335
00336
00337
00338 #ifdef TRACE_KSP
00339 double init_time = MPI_Wtime();
00340 #endif
00341
00342 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);
00343
00344 #ifdef TRACE_KSP
00345 block_diag_context->mScatterTime += MPI_Wtime() - init_time;
00346 #endif
00347
00348
00349
00350
00351
00352
00353
00354
00355 #ifdef TRACE_KSP
00356 init_time = MPI_Wtime();
00357 #endif
00358
00359 PCApply(block_diag_context->PC_amg_A11, block_diag_context->x1_subvector, block_diag_context->z);
00360 #ifdef TRACE_KSP
00361 block_diag_context->mA1PreconditionerTime += MPI_Wtime() - init_time;
00362 #endif
00363
00364
00365 #ifdef TRACE_KSP
00366 init_time = MPI_Wtime();
00367 #endif
00368 MatMult(block_diag_context->B_matrix_subblock,block_diag_context->z,block_diag_context->temp);
00369 double minus_one = -1.0;
00370 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00371 VecAYPX(&minus_one, block_diag_context->x2_subvector, block_diag_context->temp);
00372 #else
00373 VecAYPX(block_diag_context->temp, minus_one, block_diag_context->x2_subvector);
00374 #endif
00375 #ifdef TRACE_KSP
00376 block_diag_context->mExtraLAOperations += MPI_Wtime() - init_time;
00377 #endif
00378
00379
00380 #ifdef TRACE_KSP
00381 init_time = MPI_Wtime();
00382 #endif
00383 PCApply(block_diag_context->PC_amg_A22, block_diag_context->temp, block_diag_context->y2_subvector);
00384 #ifdef TRACE_KSP
00385 block_diag_context->mA2PreconditionerTime += MPI_Wtime() - init_time;
00386 #endif
00387
00388
00389
00390 #ifdef TRACE_KSP
00391 init_time = MPI_Wtime();
00392 #endif
00393 MatMult(block_diag_context->B_matrix_subblock,block_diag_context->y2_subvector,block_diag_context->temp);
00394 #ifdef TRACE_KSP
00395 block_diag_context->mExtraLAOperations += MPI_Wtime() - init_time;
00396 #endif
00397 #ifdef TRACE_KSP
00398 init_time = MPI_Wtime();
00399 #endif
00400 PCApply(block_diag_context->PC_amg_A11, block_diag_context->temp, block_diag_context->y1_subvector);
00401 #ifdef TRACE_KSP
00402 block_diag_context->mA1PreconditionerTime += MPI_Wtime() - init_time;
00403 #endif
00404
00405 #ifdef TRACE_KSP
00406 init_time = MPI_Wtime();
00407 #endif
00408 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00409 VecAYPX(&minus_one, block_diag_context->z, block_diag_context->y1_subvector);
00410 #else
00411 VecAYPX(block_diag_context->y1_subvector, minus_one, block_diag_context->z);
00412 #endif
00413 #ifdef TRACE_KSP
00414 block_diag_context->mExtraLAOperations += MPI_Wtime() - init_time;
00415 #endif
00416
00417
00418
00419
00420
00421 #ifdef TRACE_KSP
00422 init_time = MPI_Wtime();
00423 #endif
00424
00425 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);
00426
00427 #ifdef TRACE_KSP
00428 block_diag_context->mGatherTime += MPI_Wtime() - init_time;
00429 #endif
00430
00431 return 0;
00432 }