PCLDUFactorisation.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2011
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 
00029 #include <iostream>
00030 
00031 #include "PetscVecTools.hpp" // Includes Ublas so must come first
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     // Odd number of rows: impossible in Bidomain.
00095     // Odd number of local rows: impossible if V_m and phi_e for each node are stored in the same processor.
00096     if ((num_rows%2 != 0) || (num_local_rows%2 != 0))
00097     {
00098         TERMINATE("Wrong matrix parallel layout detected in PCLDUFactorisation.");
00099     }
00100 
00101     // Allocate memory
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     // Create scatter contexts
00112     {
00113         // Needed by VecScatterCreate in order to find out parallel layout.
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     // Get matrix sublock A11
00122     {
00123         // Work out local row range for subblock A11 (same as x1 or y1)
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     // Get matrix sublock A22
00146     {
00147         // Work out local row range for subblock A22 (same as x2 or y2)
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     // Get matrix sublock B (the upper triangular one)
00171     {
00172         // Work out local row range for subblock B (same as A11, x1 or y1)
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      * Experimental (#1082): in PP removing the mass matrix from the A22 block seems to work better.
00198      *                       This is equivalent to do A22 = A22 + B in this implementation.
00199      */
00200 // #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00201 //     PetscScalar petsc_one = 1.0;
00202 //     MatAXPY(&petsc_one, mPCContext.B_matrix_subblock, mPCContext.A22_matrix_subblock, DIFFERENT_NONZERO_PATTERN);
00203 // #else
00204 //     MatAXPY(mPCContext.A22_matrix_subblock, 1.0, mPCContext.B_matrix_subblock, DIFFERENT_NONZERO_PATTERN);
00205 // #endif
00206 
00207 //     // Shift the block
00208 // #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00209 //     PetscScalar shift = -1e-8;
00210 //     MatShift(&shift, mPCContext.A22_matrix_subblock);
00211 // #else
00212 //     PetscScalar shift = -1e-8;
00213 //     MatShift(mPCContext.A22_matrix_subblock, shift);
00214 // #endif
00215 
00216     PCSetType(mPetscPCObject, PCSHELL);
00217 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00218     // Register PC context and call-back function
00219     PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply, (void*) &mPCContext);
00220 #else
00221     // Register PC context so it gets passed to PCBlockDiagonalApply
00222     PCShellSetContext(mPetscPCObject, &mPCContext);
00223     // Register call-back function
00224     PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply);
00225 #endif
00226 }
00227 
00228 void PCLDUFactorisation::PCLDUFactorisationSetUp()
00229 {
00230     // These options will get read by PCSetFromOptions
00231 //     PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00232 //     PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00233 //     PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00234 
00235     /*
00236      * Set up preconditioner for block A11
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     // Choose between the two following blocks in order to approximate inv(A11) with one AMG cycle
00242     // or with an CG solve with high tolerance
00244 //    PCSetType(mPCContext.PC_amg_A11, PCBJACOBI);
00245 
00246     PCSetType(mPCContext.PC_amg_A11, PCHYPRE);
00247     //PCHYPRESetType(mPCContext.PC_amg_A11, "euclid");
00248     PetscOptionsSetValue("-pc_hypre_type", "euclid");
00249     PetscOptionsSetValue("-pc_hypre_euclid_levels", "0");
00250 
00251 //     PCSetType(mPCContext.PC_amg_A11, PCHYPRE);
00252 //     PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00253 //     PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00254 //     PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00255 //     PetscOptionsSetValue("-pc_hypre_boomeramg_coarsen_type", "HMIS");
00256 
00258 //    PCSetType(mPCContext.PC_amg_A11, PCKSP);
00259 //    KSP ksp1;
00260 //    PCKSPGetKSP(mPCContext.PC_amg_A11,&ksp1);
00261 //    KSPSetType(ksp1, KSPCG);
00262 //    KSPSetTolerances(ksp1, 0.1, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00263 //
00264 //    PC prec1;
00265 //    KSPGetPC(ksp1, &prec1);
00266 //    PCSetType(prec1, PCBJACOBI);
00267 //    PCSetFromOptions(prec1);
00268 //    PCSetOperators(prec1, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock, DIFFERENT_NONZERO_PATTERN);//   SAME_PRECONDITIONER);
00269 //    PCSetUp(prec1);
00270 //
00271 //    KSPSetFromOptions(ksp1);
00272 //    KSPSetUp(ksp1);
00274 
00275     PCSetFromOptions(mPCContext.PC_amg_A11);
00276     PCSetUp(mPCContext.PC_amg_A11);
00277 
00278     /*
00279      * Set up amg preconditioner for block A22
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     // Choose between the two following blocks in order to approximate inv(A11) with one AMG cycle
00285     // or with an CG solve with high tolerance
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     //    PetscOptionsSetValue("-pc_hypre_boomeramg_interp_type","ext+i");
00293 
00295 //    PCSetType(mPCContext.PC_amg_A22, PCKSP);
00296 //    KSP ksp2;
00297 //    PCKSPGetKSP(mPCContext.PC_amg_A22,&ksp2);
00298 //    KSPSetType(ksp2, KSPCG);
00299 //    KSPSetTolerances(ksp2, 0.1, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00300 //
00301 //    PC prec2;
00302 //    KSPGetPC(ksp2, &prec2);
00303 //    PCSetType(prec2, PCBJACOBI);
00304 //    PCSetFromOptions(prec2);
00305 //    PCSetOperators(prec2, mPCContext.A22_matrix_subblock, mPCContext.A22_matrix_subblock, DIFFERENT_NONZERO_PATTERN);//   SAME_PRECONDITIONER);
00306 //    PCSetUp(prec2);
00307 //
00308 //    KSPSetFromOptions(ksp2);
00309 //    KSPSetUp(ksp2);
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     // Cast the pointer to a PC context to our defined type
00328     PCLDUFactorisation::PCLDUFactorisationContext* block_diag_context = (PCLDUFactorisation::PCLDUFactorisationContext*) pc_context;
00329     assert(block_diag_context!=NULL);
00330 
00331     /*
00332      * Split vector x into two. x = [x1 x2]'
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      * Apply preconditioner: [y1 y2]' = inv(P)[x1 x2]'
00346      *
00347      *    z  = inv(A11)*x1
00348      *    y2 = inv(A22)*(x2 - B*z)
00349      *    y1 = z - inv(A11)(B*y2)
00350      */
00351 #ifdef TRACE_KSP
00352     init_time = MPI_Wtime();
00353 #endif
00354     //z  = inv(A11)*x1
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     //y2 = inv(A22)*(x2 - B*z)
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); //temp = B*z
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); // temp <-- x2 - temp
00368 #else
00369     VecAYPX(block_diag_context->temp, minus_one, block_diag_context->x2_subvector); // temp <-- x2 - temp
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); // y2 = inv(A22)*temp
00380 #ifdef TRACE_KSP
00381     block_diag_context->mA2PreconditionerTime += MPI_Wtime() - init_time;
00382 #endif
00383 
00384     // y1 = z - inv(A11)(B*y2)
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); //temp = B*y2
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); // y1 = inv(A11)*temp
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); // y1 <-- z - y1
00405 #else
00406     VecAYPX(block_diag_context->y1_subvector, minus_one, block_diag_context->z); // y1 <-- z - y1
00407 #endif
00408 #ifdef TRACE_KSP
00409     block_diag_context->mExtraLAOperations += MPI_Wtime() - init_time;
00410 #endif
00411 
00412     /*
00413      * Gather vectors y1 and y2. y = [y1 y2]'
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 }
Generated on Thu Dec 22 13:00:07 2011 for Chaste by  doxygen 1.6.3