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 
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     // Register PC context and call-back function
00072     PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply, (void*) &mPCContext);
00073 #else
00074     // Register PC context so it gets passed to PCBlockDiagonalApply
00075     PCShellSetContext(mPetscPCObject, &mPCContext);
00076     // Register call-back function
00077     PCShellSetApply(mPetscPCObject, PCLDUFactorisationApply);
00078 #endif
00079 
00080     // Get matrix sublock A11
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     // Get matrix sublock A22
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     // Get matrix sublock B (the upper triangular one)
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     // Allocate memory
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     // These options will get read by PCSetFromOptions
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      *  Set up preconditioner for block A11
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);//   SAME_PRECONDITIONER);
00131 
00132     // Choose between the two following blocks in order to approximate inv(A11) with one AMG cycle
00133     // or with an CG solve with high tolerance
00135     PCSetType(mPCContext.PC_amg_A11, PCHYPRE);
00137 //    PCSetType(mPCContext.PC_amg_A11, PCKSP);
00138 //    KSP ksp1;
00139 //    PCKSPGetKSP(mPCContext.PC_amg_A11,&ksp1);
00140 //    KSPSetType(ksp1, KSPCG);
00141 //    KSPSetTolerances(ksp1, 0.1, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00142 //
00143 //    PC prec1;
00144 //    KSPGetPC(ksp1, &prec1);
00145 //    PCSetType(prec1, PCBJACOBI);
00146 //    PCSetFromOptions(prec1);
00147 //    PCSetOperators(prec1, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock, DIFFERENT_NONZERO_PATTERN);//   SAME_PRECONDITIONER);
00148 //    PCSetUp(prec1);
00149 //
00150 //    KSPSetFromOptions(ksp1);
00151 //    KSPSetUp(ksp1);
00153 
00154     PCSetFromOptions(mPCContext.PC_amg_A11);
00155     PCSetUp(mPCContext.PC_amg_A11);
00156 
00157 
00158     /*
00159      *  Set up amg preconditioner for block A22
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);//   SAME_PRECONDITIONER);
00163 
00164     // Choose between the two following blocks in order to approximate inv(A11) with one AMG cycle
00165     // or with an CG solve with high tolerance
00167     PCSetType(mPCContext.PC_amg_A22, PCHYPRE);
00169 //    PCSetType(mPCContext.PC_amg_A22, PCKSP);
00170 //    KSP ksp2;
00171 //    PCKSPGetKSP(mPCContext.PC_amg_A22,&ksp2);
00172 //    KSPSetType(ksp2, KSPCG);
00173 //    KSPSetTolerances(ksp2, 0.1, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00174 //
00175 //    PC prec2;
00176 //    KSPGetPC(ksp2, &prec2);
00177 //    PCSetType(prec2, PCBJACOBI);
00178 //    PCSetFromOptions(prec2);
00179 //    PCSetOperators(prec2, mPCContext.A22_matrix_subblock, mPCContext.A22_matrix_subblock, DIFFERENT_NONZERO_PATTERN);//   SAME_PRECONDITIONER);
00180 //    PCSetUp(prec2);
00181 //
00182 //    KSPSetFromOptions(ksp2);
00183 //    KSPSetUp(ksp2);
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     // Cast the pointer to a PC context to our defined type
00195     PCLDUFactorisation::PCLDUFactorisationContext* block_diag_context = (PCLDUFactorisation::PCLDUFactorisationContext*) pc_context;
00196     assert(block_diag_context!=NULL);
00197 
00198     /*
00199      *  Split vector x into two. x = [x1 x2]'
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 //PETSc-3.x.x or PETSc-2.3.3
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 //PETSc-3.x.x or PETSc-2.3.3
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      *  Apply preconditioner: [y1 y2]' = inv(P)[x1 x2]' o
00237      *
00238      *     z  = inv(A11)*x1
00239      *     y2 = inv(A22)*(x2 - B*z)
00240      *     y1 = z - inv(A11)(B*y2)
00241      */
00242     //z  = inv(A11)*x1
00243     PCApply(block_diag_context->PC_amg_A11, block_diag_context->x1_subvector, block_diag_context->z);
00244 
00245     //y2 = inv(A22)*(x2 - B*z)
00246     MatMult(block_diag_context->B_matrix_subblock,block_diag_context->z,block_diag_context->temp); //temp = B*z
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); // temp <-- x2 - temp
00250 #else
00251     VecAYPX(block_diag_context->temp, minus_one, block_diag_context->x2_subvector); // temp <-- x2 - temp
00252 #endif
00253     PCApply(block_diag_context->PC_amg_A22, block_diag_context->temp, block_diag_context->y2_subvector); // y2 = inv(A22)*temp
00254 
00255     //y1 = z - inv(A11)(B*y2)
00256     MatMult(block_diag_context->B_matrix_subblock,block_diag_context->y2_subvector,block_diag_context->temp); //temp = B*y2
00257     PCApply(block_diag_context->PC_amg_A11, block_diag_context->temp, block_diag_context->y1_subvector); // y1 = inv(A11)*temp
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); // y1 <-- z - y1
00260 #else
00261     VecAYPX(block_diag_context->y1_subvector, minus_one, block_diag_context->z); // y1 <-- z - y1
00262 #endif
00263 
00264 
00265     /*
00266      *  Gather vectors y1 and y2. y = [y1 y2]'
00267      */
00268 //PETSc-3.x.x or PETSc-2.3.3
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 //PETSc-3.x.x or PETSc-2.3.3
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      *  Clean up
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 }

Generated by  doxygen 1.6.2