PCLDUFactorisation.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2010
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 "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     VecDestroy(mPCContext.x2_subvector);
00050     VecDestroy(mPCContext.y2_subvector);
00051     VecDestroy(mPCContext.z);
00052     VecDestroy(mPCContext.temp);
00053     
00054     VecScatterDestroy(mPCContext.A11_scatter_ctx);
00055     VecScatterDestroy(mPCContext.A22_scatter_ctx);        
00056 }
00057 
00058 void PCLDUFactorisation::PCLDUFactorisationCreate(KSP& rKspObject)
00059 {
00060     KSPGetPC(rKspObject, &mPetscPCObject);
00061 
00062     Mat system_matrix, dummy;
00063     MatStructure flag;
00064     KSPGetOperators(rKspObject, &system_matrix, &dummy, &flag);
00065 
00066     PetscInt num_rows, num_columns;
00067     MatGetSize(system_matrix, &num_rows, &num_columns);
00068 
00069     PetscInt num_local_rows, num_local_columns;
00070     MatGetLocalSize(system_matrix, &num_local_rows, &num_local_columns);
00071 
00072     // odd number of rows: impossible in Bidomain.
00073     // odd number of local rows: impossible if V_m and phi_e for each node are stored in the same processor.
00074     if ((num_rows%2 != 0) || (num_local_rows%2 != 0))
00075     {
00076         TERMINATE("Wrong matrix parallel layout detected in PCLDUFactorisation.");
00077     }
00078 
00079     // Allocate memory     
00080     unsigned subvector_num_rows = num_rows/2;    
00081     unsigned subvector_local_rows = num_local_rows/2;
00082     mPCContext.x1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00083     mPCContext.x2_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00084     mPCContext.y1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00085     mPCContext.y2_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00086     mPCContext.z = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00087     mPCContext.temp = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00088 
00089     // Create scatter contexts
00090     {
00091         // Needed by VecScatterCreate in order to find out parallel layout.
00092         Vec dummy_vec = PetscTools::CreateVec(num_rows, num_local_rows);
00093 
00094         IS A11_rows, A22_rows;        
00095         ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 2, &A11_rows);
00096         ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 1, 2, &A22_rows);   
00097         
00098         IS all_vector;    
00099         ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 1, &all_vector);
00100         
00101         VecScatterCreate(dummy_vec, A11_rows, mPCContext.x1_subvector, all_vector, &mPCContext.A11_scatter_ctx);    
00102         VecScatterCreate(dummy_vec, A22_rows, mPCContext.x2_subvector, all_vector, &mPCContext.A22_scatter_ctx);
00103         
00104         ISDestroy(A11_rows);
00105         ISDestroy(A22_rows);
00106         ISDestroy(all_vector);
00107         
00108         VecDestroy(dummy_vec);
00109     }
00110 
00111 
00112     // Get matrix sublock A11        
00113     {
00114         // Work out local row range for subblock A11 (same as x1 or y1) 
00115         PetscInt low, high, global_size;
00116         VecGetOwnershipRange(mPCContext.x1_subvector, &low, &high);
00117         VecGetSize(mPCContext.x1_subvector, &global_size);        
00118         assert(global_size == num_rows/2);       
00119         
00120         IS A11_local_rows;
00121         IS A11_columns;
00122         ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &A11_local_rows);
00123         ISCreateStride(PETSC_COMM_WORLD, global_size, 0, 2, &A11_columns);
00124     
00125 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00126         MatGetSubMatrix(system_matrix, A11_local_rows, A11_columns,
00127             MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00128 #else
00129         MatGetSubMatrix(system_matrix, A11_local_rows, A11_columns, PETSC_DECIDE,
00130             MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00131 #endif    
00132         ISDestroy(A11_local_rows);
00133         ISDestroy(A11_columns);
00134     }
00135 
00136     // Get matrix sublock A22
00137     {
00138         // Work out local row range for subblock A22 (same as x2 or y2) 
00139         PetscInt low, high, global_size;
00140         VecGetOwnershipRange(mPCContext.x2_subvector, &low, &high);
00141         VecGetSize(mPCContext.x2_subvector, &global_size);        
00142         assert(global_size == num_rows/2);       
00143         
00144         IS A22_local_rows;
00145         IS A22_columns;
00146         ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low+1, 2, &A22_local_rows);
00147         ISCreateStride(PETSC_COMM_WORLD, global_size, 1, 2, &A22_columns);
00148 
00149 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00150         MatGetSubMatrix(system_matrix, A22_local_rows, A22_columns,
00151             MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00152 #else
00153         MatGetSubMatrix(system_matrix, A22_local_rows, A22_columns, PETSC_DECIDE, 
00154             MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00155 #endif
00156     
00157         ISDestroy(A22_local_rows);
00158         ISDestroy(A22_columns);
00159     }
00160 
00161     // Get matrix sublock B (the upper triangular one)
00162     {
00163         // Work out local row range for subblock B (same as A11, x1 or y1) 
00164         PetscInt low, high, global_size;
00165         VecGetOwnershipRange(mPCContext.x1_subvector, &low, &high);
00166         VecGetSize(mPCContext.x1_subvector, &global_size);        
00167         assert(global_size == num_rows/2);        
00168         
00169         IS B_local_rows;
00170         IS B_columns;
00171         ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &B_local_rows);
00172         ISCreateStride(PETSC_COMM_WORLD, global_size, 1, 2, &B_columns);
00173     
00174 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00175         MatGetSubMatrix(system_matrix, B_local_rows, B_columns,
00176             MAT_INITIAL_MATRIX, &mPCContext.B_matrix_subblock);        
00177 #else
00178         MatGetSubMatrix(system_matrix, B_local_rows, B_columns, PETSC_DECIDE, 
00179             MAT_INITIAL_MATRIX, &mPCContext.B_matrix_subblock);
00180 #endif
00181     
00182         ISDestroy(B_local_rows);
00183         ISDestroy(B_columns);
00184     }
00185     
00186     /*
00187      * Experimental (#1082): in PP removing the mass matrix from the A22 block seems to work better.
00188      *                       This is equivalent to do A22 = A22 + B in this implementation. 
00189      */
00190 // #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00191 //     PetscScalar petsc_one = 1.0;
00192 //     MatAXPY(&petsc_one, mPCContext.B_matrix_subblock, mPCContext.A22_matrix_subblock, DIFFERENT_NONZERO_PATTERN);
00193 // #else
00194 //     MatAXPY(mPCContext.A22_matrix_subblock, 1.0, mPCContext.B_matrix_subblock, DIFFERENT_NONZERO_PATTERN);
00195 // #endif    
00196     
00197 //     // Shift the block
00198 // #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00199 //     PetscScalar shift = -1e-8;
00200 //     MatShift(&shift, mPCContext.A22_matrix_subblock);
00201 // #else
00202 //     PetscScalar shift = -1e-8;
00203 //     MatShift(mPCContext.A22_matrix_subblock, shift);
00204 // #endif    
00205 
00206 
00207 
00208     PCSetType(mPetscPCObject, PCSHELL);
00209 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00210     // Register PC context and call-back function
00211     PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply, (void*) &mPCContext);
00212 #else
00213     // Register PC context so it gets passed to PCBlockDiagonalApply
00214     PCShellSetContext(mPetscPCObject, &mPCContext);
00215     // Register call-back function
00216     PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply);
00217 #endif
00218 
00219 }
00220 
00221 void PCLDUFactorisation::PCLDUFactorisationSetUp()
00222 {
00223     // These options will get read by PCSetFromOptions
00224     PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00225     PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00226     PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00227 
00228     /*
00229      *  Set up preconditioner for block A11
00230      */
00231     PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A11));
00232     PCSetOperators(mPCContext.PC_amg_A11, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock, DIFFERENT_NONZERO_PATTERN);//   SAME_PRECONDITIONER);
00233 
00234     // Choose between the two following blocks in order to approximate inv(A11) with one AMG cycle
00235     // or with an CG solve with high tolerance
00237     PCSetType(mPCContext.PC_amg_A11, PCBJACOBI);
00238     //PCSetType(mPCContext.PC_amg_A11, PCHYPRE);
00239     //PCHYPRESetType(mPCContext.PC_amg_A11, "euclid");
00240 
00242 //    PCSetType(mPCContext.PC_amg_A11, PCKSP);
00243 //    KSP ksp1;
00244 //    PCKSPGetKSP(mPCContext.PC_amg_A11,&ksp1);
00245 //    KSPSetType(ksp1, KSPCG);
00246 //    KSPSetTolerances(ksp1, 0.1, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00247 //
00248 //    PC prec1;
00249 //    KSPGetPC(ksp1, &prec1);
00250 //    PCSetType(prec1, PCBJACOBI);
00251 //    PCSetFromOptions(prec1);
00252 //    PCSetOperators(prec1, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock, DIFFERENT_NONZERO_PATTERN);//   SAME_PRECONDITIONER);
00253 //    PCSetUp(prec1);
00254 //
00255 //    KSPSetFromOptions(ksp1);
00256 //    KSPSetUp(ksp1);
00258 
00259     PCSetFromOptions(mPCContext.PC_amg_A11);
00260     PCSetUp(mPCContext.PC_amg_A11);
00261 
00262 
00263     /*
00264      *  Set up amg preconditioner for block A22
00265      */
00266     PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A22));
00267     PCSetOperators(mPCContext.PC_amg_A22, mPCContext.A22_matrix_subblock, mPCContext.A22_matrix_subblock, DIFFERENT_NONZERO_PATTERN);//   SAME_PRECONDITIONER);
00268 
00269     // Choose between the two following blocks in order to approximate inv(A11) with one AMG cycle
00270     // or with an CG solve with high tolerance
00272     PCSetType(mPCContext.PC_amg_A22, PCHYPRE);
00274 //    PCSetType(mPCContext.PC_amg_A22, PCKSP);
00275 //    KSP ksp2;
00276 //    PCKSPGetKSP(mPCContext.PC_amg_A22,&ksp2);
00277 //    KSPSetType(ksp2, KSPCG);
00278 //    KSPSetTolerances(ksp2, 0.1, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00279 //
00280 //    PC prec2;
00281 //    KSPGetPC(ksp2, &prec2);
00282 //    PCSetType(prec2, PCBJACOBI);
00283 //    PCSetFromOptions(prec2);
00284 //    PCSetOperators(prec2, mPCContext.A22_matrix_subblock, mPCContext.A22_matrix_subblock, DIFFERENT_NONZERO_PATTERN);//   SAME_PRECONDITIONER);
00285 //    PCSetUp(prec2);
00286 //
00287 //    KSPSetFromOptions(ksp2);
00288 //    KSPSetUp(ksp2);
00290 
00291     PCSetFromOptions(mPCContext.PC_amg_A22);
00292     PCSetUp(mPCContext.PC_amg_A22);
00293 }
00294 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00295 PetscErrorCode PCLDUFactorisationApply(PC pc_object, Vec x, Vec y)
00296 {
00297   void* pc_context;
00298 
00299   PCShellGetContext(pc_object, &pc_context);   
00300 #else
00301 PetscErrorCode PCLDUFactorisationApply(void* pc_context, Vec x, Vec y)
00302 {
00303 #endif
00305 
00306     // Cast the pointer to a PC context to our defined type
00307     PCLDUFactorisation::PCLDUFactorisationContext* block_diag_context = (PCLDUFactorisation::PCLDUFactorisationContext*) pc_context;
00308     assert(block_diag_context!=NULL);
00309 
00310     /*
00311      *  Split vector x into two. x = [x1 x2]'
00312      */
00313 //PETSc-3.x.x or PETSc-2.3.3
00314 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00315     VecScatterBegin(block_diag_context->A11_scatter_ctx, x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD);
00316     VecScatterEnd(block_diag_context->A11_scatter_ctx, x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD);
00317 #else
00318     VecScatterBegin(x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A11_scatter_ctx);
00319     VecScatterEnd(x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A11_scatter_ctx);
00320 #endif
00321 
00322 //PETSc-3.x.x or PETSc-2.3.3
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->A22_scatter_ctx, x, block_diag_context->x2_subvector, INSERT_VALUES, SCATTER_FORWARD);
00325     VecScatterEnd(block_diag_context->A22_scatter_ctx, x, block_diag_context->x2_subvector, INSERT_VALUES, SCATTER_FORWARD);
00326 #else
00327     VecScatterBegin(x, block_diag_context->x2_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A22_scatter_ctx);
00328     VecScatterEnd(x, block_diag_context->x2_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A22_scatter_ctx);
00329 #endif
00330 
00331 
00332     /*
00333      *  Apply preconditioner: [y1 y2]' = inv(P)[x1 x2]' 
00334      *
00335      *     z  = inv(A11)*x1
00336      *     y2 = inv(A22)*(x2 - B*z)
00337      *     y1 = z - inv(A11)(B*y2)
00338      */
00339     //z  = inv(A11)*x1
00340     PCApply(block_diag_context->PC_amg_A11, block_diag_context->x1_subvector, block_diag_context->z);
00341 
00342     //y2 = inv(A22)*(x2 - B*z)
00343     MatMult(block_diag_context->B_matrix_subblock,block_diag_context->z,block_diag_context->temp); //temp = B*z
00344     double minus_one = -1.0;
00345 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00346     VecAYPX(&minus_one, block_diag_context->x2_subvector, block_diag_context->temp); // temp <-- x2 - temp
00347 #else
00348     VecAYPX(block_diag_context->temp, minus_one, block_diag_context->x2_subvector); // temp <-- x2 - temp
00349 #endif
00350     PCApply(block_diag_context->PC_amg_A22, block_diag_context->temp, block_diag_context->y2_subvector); // y2 = inv(A22)*temp
00351 
00352     //y1 = z - inv(A11)(B*y2)
00353     MatMult(block_diag_context->B_matrix_subblock,block_diag_context->y2_subvector,block_diag_context->temp); //temp = B*y2
00354     PCApply(block_diag_context->PC_amg_A11, block_diag_context->temp, block_diag_context->y1_subvector); // y1 = inv(A11)*temp
00355 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00356     VecAYPX(&minus_one, block_diag_context->z, block_diag_context->y1_subvector); // y1 <-- z - y1
00357 #else
00358     VecAYPX(block_diag_context->y1_subvector, minus_one, block_diag_context->z); // y1 <-- z - y1
00359 #endif
00360 
00361 
00362     /*
00363      *  Gather vectors y1 and y2. y = [y1 y2]'
00364      */
00365 //PETSc-3.x.x or PETSc-2.3.3
00366 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00367     VecScatterBegin(block_diag_context->A11_scatter_ctx, block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00368     VecScatterEnd(block_diag_context->A11_scatter_ctx, block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00369 #else
00370     VecScatterBegin(block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A11_scatter_ctx);
00371     VecScatterEnd(block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A11_scatter_ctx);
00372 #endif
00373 
00374 //PETSc-3.x.x or PETSc-2.3.3
00375 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00376     VecScatterBegin(block_diag_context->A22_scatter_ctx, block_diag_context->y2_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00377     VecScatterEnd(block_diag_context->A22_scatter_ctx, block_diag_context->y2_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00378 #else
00379     VecScatterBegin(block_diag_context->y2_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A22_scatter_ctx);
00380     VecScatterEnd(block_diag_context->y2_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A22_scatter_ctx);
00381 #endif
00382 
00383     return 0;
00384 }

Generated on Mon Nov 1 12:35:17 2010 for Chaste by  doxygen 1.5.5