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 PetscInt low, high, global_size;
00125 VecGetOwnershipRange(mPCContext.x1_subvector, &low, &high);
00126 VecGetSize(mPCContext.x1_subvector, &global_size);
00127 assert(global_size == num_rows/2);
00128
00129 IS A11_local_rows;
00130 IS A11_columns;
00131 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &A11_local_rows);
00132 ISCreateStride(PETSC_COMM_WORLD, global_size, 0, 2, &A11_columns);
00133
00134 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00135 MatGetSubMatrix(system_matrix, A11_local_rows, A11_local_rows,
00136 MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00137 #else
00138 MatGetSubMatrix(system_matrix, A11_local_rows, A11_columns, PETSC_DECIDE,
00139 MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00140 #endif
00141 ISDestroy(A11_local_rows);
00142 ISDestroy(A11_columns);
00143 }
00144
00145
00146 {
00147
00148 PetscInt low, high, global_size;
00149 VecGetOwnershipRange(mPCContext.x2_subvector, &low, &high);
00150 VecGetSize(mPCContext.x2_subvector, &global_size);
00151 assert(global_size == num_rows/2);
00152
00153 IS A22_local_rows;
00154 IS A22_columns;
00155 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low+1, 2, &A22_local_rows);
00156 ISCreateStride(PETSC_COMM_WORLD, global_size, 1, 2, &A22_columns);
00157
00158 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00159 MatGetSubMatrix(system_matrix, A22_local_rows, A22_local_rows,
00160 MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00161 #else
00162 MatGetSubMatrix(system_matrix, A22_local_rows, A22_columns, PETSC_DECIDE,
00163 MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00164 #endif
00165
00166 ISDestroy(A22_local_rows);
00167 ISDestroy(A22_columns);
00168 }
00169
00170
00171 {
00172
00173 PetscInt low, high, global_size;
00174 VecGetOwnershipRange(mPCContext.x1_subvector, &low, &high);
00175 VecGetSize(mPCContext.x1_subvector, &global_size);
00176 assert(global_size == num_rows/2);
00177
00178 IS B_local_rows;
00179 IS B_columns;
00180 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &B_local_rows);
00181
00182 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00183 ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low+1, 2, &B_columns);
00184 MatGetSubMatrix(system_matrix, B_local_rows, B_columns,
00185 MAT_INITIAL_MATRIX, &mPCContext.B_matrix_subblock);
00186 #else
00187 ISCreateStride(PETSC_COMM_WORLD, global_size, 1, 2, &B_columns);
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 PCSetType(mPetscPCObject, PCSHELL);
00217 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00218
00219 PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply, (void*) &mPCContext);
00220 #else
00221
00222 PCShellSetContext(mPetscPCObject, &mPCContext);
00223
00224 PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply);
00225 #endif
00226 }
00227
00228 void PCLDUFactorisation::PCLDUFactorisationSetUp()
00229 {
00230
00231
00232
00233
00234
00235
00236
00237
00238 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A11));
00239 PCSetOperators(mPCContext.PC_amg_A11, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock, SAME_PRECONDITIONER);
00240
00241
00242
00244
00245
00246 PCSetType(mPCContext.PC_amg_A11, PCHYPRE);
00247
00248 PetscOptionsSetValue("-pc_hypre_type", "euclid");
00249 PetscOptionsSetValue("-pc_hypre_euclid_levels", "0");
00250
00251
00252
00253
00254
00255
00256
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272
00274
00275 PCSetFromOptions(mPCContext.PC_amg_A11);
00276 PCSetUp(mPCContext.PC_amg_A11);
00277
00278
00279
00280
00281 PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A22));
00282 PCSetOperators(mPCContext.PC_amg_A22, mPCContext.A22_matrix_subblock, mPCContext.A22_matrix_subblock, SAME_PRECONDITIONER);
00283
00284
00285
00287 PCSetType(mPCContext.PC_amg_A22, PCHYPRE);
00288 PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00289 PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00290 PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00291 PetscOptionsSetValue("-pc_hypre_boomeramg_coarsen_type", "HMIS");
00292
00293
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00311
00312 PCSetFromOptions(mPCContext.PC_amg_A22);
00313 PCSetUp(mPCContext.PC_amg_A22);
00314 }
00315 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00316 PetscErrorCode PCLDUFactorisationApply(PC pc_object, Vec x, Vec y)
00317 {
00318 void* pc_context;
00319
00320 PCShellGetContext(pc_object, &pc_context);
00321 #else
00322 PetscErrorCode PCLDUFactorisationApply(void* pc_context, Vec x, Vec y)
00323 {
00324 #endif
00325
00326
00327
00328 PCLDUFactorisation::PCLDUFactorisationContext* block_diag_context = (PCLDUFactorisation::PCLDUFactorisationContext*) pc_context;
00329 assert(block_diag_context!=NULL);
00330
00331
00332
00333
00334 #ifdef TRACE_KSP
00335 double init_time = MPI_Wtime();
00336 #endif
00337
00338 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);
00339
00340 #ifdef TRACE_KSP
00341 block_diag_context->mScatterTime += MPI_Wtime() - init_time;
00342 #endif
00343
00344
00345
00346
00347
00348
00349
00350
00351 #ifdef TRACE_KSP
00352 init_time = MPI_Wtime();
00353 #endif
00354
00355 PCApply(block_diag_context->PC_amg_A11, block_diag_context->x1_subvector, block_diag_context->z);
00356 #ifdef TRACE_KSP
00357 block_diag_context->mA1PreconditionerTime += MPI_Wtime() - init_time;
00358 #endif
00359
00360
00361 #ifdef TRACE_KSP
00362 init_time = MPI_Wtime();
00363 #endif
00364 MatMult(block_diag_context->B_matrix_subblock,block_diag_context->z,block_diag_context->temp);
00365 double minus_one = -1.0;
00366 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00367 VecAYPX(&minus_one, block_diag_context->x2_subvector, block_diag_context->temp);
00368 #else
00369 VecAYPX(block_diag_context->temp, minus_one, block_diag_context->x2_subvector);
00370 #endif
00371 #ifdef TRACE_KSP
00372 block_diag_context->mExtraLAOperations += MPI_Wtime() - init_time;
00373 #endif
00374
00375
00376 #ifdef TRACE_KSP
00377 init_time = MPI_Wtime();
00378 #endif
00379 PCApply(block_diag_context->PC_amg_A22, block_diag_context->temp, block_diag_context->y2_subvector);
00380 #ifdef TRACE_KSP
00381 block_diag_context->mA2PreconditionerTime += MPI_Wtime() - init_time;
00382 #endif
00383
00384
00385 #ifdef TRACE_KSP
00386 init_time = MPI_Wtime();
00387 #endif
00388 MatMult(block_diag_context->B_matrix_subblock,block_diag_context->y2_subvector,block_diag_context->temp);
00389 #ifdef TRACE_KSP
00390 block_diag_context->mExtraLAOperations += MPI_Wtime() - init_time;
00391 #endif
00392 #ifdef TRACE_KSP
00393 init_time = MPI_Wtime();
00394 #endif
00395 PCApply(block_diag_context->PC_amg_A11, block_diag_context->temp, block_diag_context->y1_subvector);
00396 #ifdef TRACE_KSP
00397 block_diag_context->mA1PreconditionerTime += MPI_Wtime() - init_time;
00398 #endif
00399
00400 #ifdef TRACE_KSP
00401 init_time = MPI_Wtime();
00402 #endif
00403 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00404 VecAYPX(&minus_one, block_diag_context->z, block_diag_context->y1_subvector);
00405 #else
00406 VecAYPX(block_diag_context->y1_subvector, minus_one, block_diag_context->z);
00407 #endif
00408 #ifdef TRACE_KSP
00409 block_diag_context->mExtraLAOperations += MPI_Wtime() - init_time;
00410 #endif
00411
00412
00413
00414
00415 #ifdef TRACE_KSP
00416 init_time = MPI_Wtime();
00417 #endif
00418
00419 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);
00420
00421 #ifdef TRACE_KSP
00422 block_diag_context->mGatherTime += MPI_Wtime() - init_time;
00423 #endif
00424
00425 return 0;
00426 }