PCTwoLevelsBlockDiagonal.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 "PCTwoLevelsBlockDiagonal.hpp"
00030 #include "Exception.hpp"
00031 
00032 #include <iostream>
00033 
00034 PCTwoLevelsBlockDiagonal::PCTwoLevelsBlockDiagonal(KSP& rKspObject, std::vector<PetscInt>& rBathNodes)
00035 {
00036     PCTwoLevelsBlockDiagonalCreate(rKspObject, rBathNodes);
00037     PCTwoLevelsBlockDiagonalSetUp();
00038 }
00039 
00040 PCTwoLevelsBlockDiagonal::~PCTwoLevelsBlockDiagonal()
00041 {
00042     MatDestroy(mPCContext.A11_matrix_subblock);
00043     MatDestroy(mPCContext.A22_B1_matrix_subblock);
00044     MatDestroy(mPCContext.A22_B2_matrix_subblock);
00045 
00046     PCDestroy(mPCContext.PC_amg_A11);
00047     PCDestroy(mPCContext.PC_amg_A22_B1);
00048     PCDestroy(mPCContext.PC_amg_A22_B2);
00049 
00050     VecDestroy(mPCContext.x1_subvector);
00051     VecDestroy(mPCContext.y1_subvector);
00052 
00053     VecDestroy(mPCContext.x21_subvector);
00054     VecDestroy(mPCContext.y21_subvector);
00055 
00056     VecDestroy(mPCContext.x22_subvector);
00057     VecDestroy(mPCContext.y22_subvector);
00058 
00059     VecScatterDestroy(mPCContext.A11_scatter_ctx);
00060     VecScatterDestroy(mPCContext.A22_B1_scatter_ctx);
00061     VecScatterDestroy(mPCContext.A22_B2_scatter_ctx);
00062 }
00063 
00064 void PCTwoLevelsBlockDiagonal::PCTwoLevelsBlockDiagonalCreate(KSP& rKspObject, std::vector<PetscInt>& rBathNodes)
00065 {
00066     KSPGetPC(rKspObject, &mPetscPCObject);
00067 
00068     Mat system_matrix, dummy;
00069     MatStructure flag;
00070     KSPGetOperators(rKspObject, &system_matrix, &dummy, &flag);
00071 
00072     PetscInt num_rows, num_columns;
00073     MatGetSize(system_matrix, &num_rows, &num_columns);
00074     assert(num_rows==num_columns);
00075 
00076     PetscInt num_local_rows, num_local_columns;
00077     MatGetLocalSize(system_matrix, &num_local_rows, &num_local_columns);
00078 
00079     // Odd number of rows: impossible in Bidomain.
00080     // Odd number of local rows: impossible if V_m and phi_e for each node are stored in the same processor.
00081     if ((num_rows%2 != 0) || (num_local_rows%2 != 0))
00082     {
00083         TERMINATE("Wrong matrix parallel layout detected in PCLDUFactorisation.");
00084     }
00085 
00086     // Allocate memory
00087     unsigned subvector_num_rows = num_rows/2;
00088     unsigned subvector_local_rows = num_local_rows/2;
00089 
00090     unsigned subvector_num_rows_tissue = subvector_num_rows - rBathNodes.size();
00091     unsigned subvector_local_rows_tissue = subvector_num_rows_tissue; 
00092 
00093     unsigned subvector_num_rows_bath = rBathNodes.size();
00094     unsigned subvector_local_rows_bath = subvector_num_rows_bath; 
00095 
00096     assert(PetscTools::IsSequential());
00097 
00098     mPCContext.x1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00099     mPCContext.x21_subvector = PetscTools::CreateVec(subvector_num_rows_tissue, subvector_local_rows_tissue);
00100     mPCContext.x22_subvector = PetscTools::CreateVec(subvector_num_rows_bath, subvector_local_rows_bath);
00101     mPCContext.y1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows);
00102     mPCContext.y21_subvector = PetscTools::CreateVec(subvector_num_rows_tissue, subvector_local_rows_tissue);
00103     mPCContext.y22_subvector = PetscTools::CreateVec(subvector_num_rows_bath, subvector_local_rows_bath);
00104 
00105     // Define IS objects that will be used throughout the method
00106     IS A11_all_rows;
00107     ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 2, &A11_all_rows);
00108 
00109     IS A22_all_rows;
00110     PetscInt A22_size;
00111     VecGetSize(mPCContext.x1_subvector, &A22_size);
00113     ISCreateStride(PETSC_COMM_WORLD, A22_size, 1, 2, &A22_all_rows);
00114 
00115     IS A22_bath_rows;
00116     PetscInt* phi_e_bath_rows = new PetscInt[rBathNodes.size()];
00117     for (unsigned index=0; index<rBathNodes.size(); index++)
00118     {
00119         phi_e_bath_rows[index] = 2*rBathNodes[index] + 1;
00120     }
00121     ISCreateGeneralWithArray(PETSC_COMM_WORLD, rBathNodes.size(), phi_e_bath_rows, &A22_bath_rows);
00122 
00123     IS A22_tissue_rows;
00124     ISDifference(A22_all_rows, A22_bath_rows, &A22_tissue_rows);
00125 
00126     // Create scatter contexts
00127     {
00128         // Needed by VecScatterCreate in order to find out parallel layout.
00129         Vec dummy_vec = PetscTools::CreateVec(num_rows, num_local_rows);
00130 
00132         IS& A11_rows=A11_all_rows;
00133         IS& A22_B1_rows=A22_tissue_rows;
00134         IS& A22_B2_rows=A22_bath_rows;
00135 
00136         IS all_vector;
00137         ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 1, &all_vector);
00138 
00139         IS tissue_vector;
00140         ISCreateStride(PETSC_COMM_WORLD, (num_rows/2)-rBathNodes.size(), 0, 1, &tissue_vector);
00141 
00142         IS bath_vector;
00143         ISCreateStride(PETSC_COMM_WORLD, rBathNodes.size(), 0, 1, &bath_vector);
00144 
00145         VecScatterCreate(dummy_vec, A11_rows, mPCContext.x1_subvector, all_vector, &mPCContext.A11_scatter_ctx);
00146         VecScatterCreate(dummy_vec, A22_B1_rows, mPCContext.x21_subvector, tissue_vector, &mPCContext.A22_B1_scatter_ctx);
00147         VecScatterCreate(dummy_vec, A22_B2_rows, mPCContext.x22_subvector, bath_vector, &mPCContext.A22_B2_scatter_ctx);
00148 
00149         ISDestroy(all_vector);
00150         ISDestroy(tissue_vector);
00151         ISDestroy(bath_vector);
00152 
00153         VecDestroy(dummy_vec);
00154     }
00155 
00156     // Get matrix sublock A11
00157     {
00158         // Work out local row range for subblock A11 (same as x1 or y1)
00159         PetscInt low, high, global_size;
00160         VecGetOwnershipRange(mPCContext.x1_subvector, &low, &high);
00161         VecGetSize(mPCContext.x1_subvector, &global_size);
00162         assert(global_size == num_rows/2);
00163 
00164         IS A11_local_rows;
00165         IS& A11_columns=A11_all_rows;
00166         ISCreateStride(PETSC_COMM_WORLD, high-low, 2*low, 2, &A11_local_rows); 
00167 
00168 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00169         MatGetSubMatrix(system_matrix, A11_local_rows, A11_columns,
00170             MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00171 #else
00172         MatGetSubMatrix(system_matrix, A11_local_rows, A11_columns, PETSC_DECIDE,
00173             MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00174 #endif
00175 
00176         ISDestroy(A11_local_rows);
00177     }
00178 
00179     // Get matrix sublock A22_B1
00180     {
00181 //        // Work out local row range for subblock A22 (same as x2 or y2)
00182 //        PetscInt low, high, global_size;
00183 //        VecGetOwnershipRange(mPCContext.x21_subvector, &low, &high);
00184 //        VecGetSize(mPCContext.x21_subvector, &global_size);
00185 //        assert(global_size == (num_rows/2) - (PetscInt) rBathNodes.size());
00186 
00187         assert(PetscTools::IsSequential());
00188         IS& A22_B1_local_rows = A22_tissue_rows; // wrong in parallel, need to give local rows
00189         IS& A22_B1_columns = A22_tissue_rows;
00190 
00191 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00192         MatGetSubMatrix(system_matrix, A22_B1_local_rows, A22_B1_columns,
00193             MAT_INITIAL_MATRIX, &mPCContext.A22_B1_matrix_subblock);
00194 #else
00195         MatGetSubMatrix(system_matrix, A22_B1_local_rows, A22_B1_columns, PETSC_DECIDE,
00196             MAT_INITIAL_MATRIX, &mPCContext.A22_B1_matrix_subblock);
00197 #endif
00198 
00199     }
00200 
00201     // Get matrix sublock A22_B2
00202     {
00203 //        // Work out local row range for subblock A22 (same as x2 or y2)
00204 //        PetscInt low, high, global_size;
00205 //        VecGetOwnershipRange(mPCContext.x21_subvector, &low, &high);
00206 //        VecGetSize(mPCContext.x21_subvector, &global_size);
00207 //        assert(global_size == (num_rows/2) - (PetscInt) rBathNodes.size());
00208 
00209         assert(PetscTools::IsSequential());
00210         IS& A22_B2_local_rows = A22_bath_rows; // wrong in parallel, need to give local rows
00211         IS& A22_B2_columns = A22_bath_rows;
00212 
00213 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00214         MatGetSubMatrix(system_matrix, A22_B2_local_rows, A22_B2_columns,
00215             MAT_INITIAL_MATRIX, &mPCContext.A22_B2_matrix_subblock);
00216 #else
00217         MatGetSubMatrix(system_matrix, A22_B2_local_rows, A22_B2_columns, PETSC_DECIDE,
00218             MAT_INITIAL_MATRIX, &mPCContext.A22_B2_matrix_subblock);
00219 #endif
00220     }
00221 
00222     ISDestroy(A11_all_rows);
00223     ISDestroy(A22_all_rows);
00224     ISDestroy(A22_bath_rows);
00225     delete[] phi_e_bath_rows;
00226     ISDestroy(A22_tissue_rows);
00227 
00228     // Register call-back function and its context
00229     PCSetType(mPetscPCObject, PCSHELL);
00230 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00231     PCShellSetApply(mPetscPCObject, PCTwoLevelsBlockDiagonalApply, (void*) &mPCContext);
00232 #else
00233     // Register PC context so it gets passed to PCTwoLevelsBlockDiagonalApply
00234     PCShellSetContext(mPetscPCObject, &mPCContext);
00235 
00236     // Register call-back function
00237     PCShellSetApply(mPetscPCObject, PCTwoLevelsBlockDiagonalApply);
00238 #endif
00239 }
00240 
00241 void PCTwoLevelsBlockDiagonal::PCTwoLevelsBlockDiagonalSetUp()
00242 {
00243     // These options will get read by PCSetFromOptions
00244     PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00245     PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00246     PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00247 
00248     // Set up amg preconditioner for block A11
00249     PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A11));
00250     PCSetType(mPCContext.PC_amg_A11, PCBJACOBI);
00251     PCSetOperators(mPCContext.PC_amg_A11, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock, DIFFERENT_NONZERO_PATTERN);//   SAME_PRECONDITIONER);
00252     PCSetFromOptions(mPCContext.PC_amg_A11);
00253     PCSetUp(mPCContext.PC_amg_A11);
00254 
00255     // Set up amg preconditioner for block A22_B1
00256     PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A22_B1));
00257     PCSetType(mPCContext.PC_amg_A22_B1, PCBJACOBI);
00258     PCSetOperators(mPCContext.PC_amg_A22_B1, mPCContext.A22_B1_matrix_subblock, mPCContext.A22_B1_matrix_subblock, DIFFERENT_NONZERO_PATTERN);//   SAME_PRECONDITIONER);
00259     PCSetFromOptions(mPCContext.PC_amg_A22_B1);
00260     PCSetUp(mPCContext.PC_amg_A22_B1);
00261 
00262     // Set up amg preconditioner for block A22_B2
00263     PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A22_B2));
00264     PCSetType(mPCContext.PC_amg_A22_B2, PCHYPRE);
00265     //PCHYPRESetType(mPCContext.PC_amg_A22_B2, "boomeramg");
00266     PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00267 
00268     PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00269     PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00270     PetscOptionsSetValue("-pc_hypre_boomeramg_coarsen_type", "HMIS");
00271 
00272     PCSetOperators(mPCContext.PC_amg_A22_B2, mPCContext.A22_B2_matrix_subblock, mPCContext.A22_B2_matrix_subblock, DIFFERENT_NONZERO_PATTERN);//   SAME_PRECONDITIONER);
00273     PCSetFromOptions(mPCContext.PC_amg_A22_B2);
00274     PCSetUp(mPCContext.PC_amg_A22_B2);
00275 }
00276 
00277 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 1) //PETSc 3.1
00278 PetscErrorCode PCTwoLevelsBlockDiagonalApply(PC pc_object, Vec x, Vec y)
00279 {
00280   void* pc_context;
00281 
00282   PCShellGetContext(pc_object, &pc_context);
00283 #else
00284 PetscErrorCode PCTwoLevelsBlockDiagonalApply(void* pc_context, Vec x, Vec y)
00285 {
00286 #endif
00287 
00288     // Cast the context pointer to PCTwoLevelsBlockDiagonalContext
00289     PCTwoLevelsBlockDiagonal::PCTwoLevelsBlockDiagonalContext* block_diag_context = (PCTwoLevelsBlockDiagonal::PCTwoLevelsBlockDiagonalContext*) pc_context;
00290     assert(block_diag_context!=NULL);
00291 
00292     /*
00293      * Scatter x = [x1 x21 x22]'
00294      */
00295 //PETSc-3.x.x or PETSc-2.3.3
00296 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00297     VecScatterBegin(block_diag_context->A11_scatter_ctx, x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD);
00298     VecScatterEnd(block_diag_context->A11_scatter_ctx, x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD);
00299 #else
00300     VecScatterBegin(x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A11_scatter_ctx);
00301     VecScatterEnd(x, block_diag_context->x1_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A11_scatter_ctx);
00302 #endif
00303 
00304 //PETSc-3.x.x or PETSc-2.3.3
00305 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00306     VecScatterBegin(block_diag_context->A22_B1_scatter_ctx, x, block_diag_context->x21_subvector, INSERT_VALUES, SCATTER_FORWARD);
00307     VecScatterEnd(block_diag_context->A22_B1_scatter_ctx, x, block_diag_context->x21_subvector, INSERT_VALUES, SCATTER_FORWARD);
00308 #else
00309     VecScatterBegin(x, block_diag_context->x21_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A22_B1_scatter_ctx);
00310     VecScatterEnd(x, block_diag_context->x21_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A22_B1_scatter_ctx);
00311 #endif
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->A22_B2_scatter_ctx, x, block_diag_context->x22_subvector, INSERT_VALUES, SCATTER_FORWARD);
00316     VecScatterEnd(block_diag_context->A22_B2_scatter_ctx, x, block_diag_context->x22_subvector, INSERT_VALUES, SCATTER_FORWARD);
00317 #else
00318     VecScatterBegin(x, block_diag_context->x22_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A22_B2_scatter_ctx);
00319     VecScatterEnd(x, block_diag_context->x22_subvector, INSERT_VALUES, SCATTER_FORWARD, block_diag_context->A22_B2_scatter_ctx);
00320 #endif
00321 
00322     /*
00323      * y1  = ILU(A11)*x1
00324      * y21 = ILU(A22)*x21
00325      * y22 = AMG(A22)*x22
00326      */
00327     PCApply(block_diag_context->PC_amg_A11, block_diag_context->x1_subvector, block_diag_context->y1_subvector);
00328     PCApply(block_diag_context->PC_amg_A22_B1, block_diag_context->x21_subvector, block_diag_context->y21_subvector);
00329     PCApply(block_diag_context->PC_amg_A22_B2, block_diag_context->x22_subvector, block_diag_context->y22_subvector);
00330 
00331     /*
00332      * Gather y = [y1 y21 y22]'
00333      */
00334 //PETSc-3.x.x or PETSc-2.3.3
00335 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00336     VecScatterBegin(block_diag_context->A11_scatter_ctx, block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00337     VecScatterEnd(block_diag_context->A11_scatter_ctx, block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00338 #else
00339     VecScatterBegin(block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A11_scatter_ctx);
00340     VecScatterEnd(block_diag_context->y1_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A11_scatter_ctx);
00341 #endif
00342 
00343 //PETSc-3.x.x or PETSc-2.3.3
00344 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00345     VecScatterBegin(block_diag_context->A22_B1_scatter_ctx, block_diag_context->y21_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00346     VecScatterEnd(block_diag_context->A22_B1_scatter_ctx, block_diag_context->y21_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00347 #else
00348     VecScatterBegin(block_diag_context->y21_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A22_B1_scatter_ctx);
00349     VecScatterEnd(block_diag_context->y21_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A22_B1_scatter_ctx);
00350 #endif
00351 
00352 //PETSc-3.x.x or PETSc-2.3.3
00353 #if ( (PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
00354     VecScatterBegin(block_diag_context->A22_B2_scatter_ctx, block_diag_context->y22_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00355     VecScatterEnd(block_diag_context->A22_B2_scatter_ctx, block_diag_context->y22_subvector, y, INSERT_VALUES, SCATTER_REVERSE);
00356 #else
00357     VecScatterBegin(block_diag_context->y22_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A22_B2_scatter_ctx);
00358     VecScatterEnd(block_diag_context->y22_subvector, y, INSERT_VALUES, SCATTER_REVERSE, block_diag_context->A22_B2_scatter_ctx);
00359 #endif
00360 
00361     return 0;
00362 }
Generated on Thu Dec 22 13:00:07 2011 for Chaste by  doxygen 1.6.3