PCBlockDiagonal.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2009
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 "PCBlockDiagonal.hpp"
00030 
00031 PCBlockDiagonal::PCBlockDiagonal(KSP& rKspObject)
00032 {    
00033     PCBlockDiagonalCreate(rKspObject);    
00034     PCBlockDiagonalSetUp();
00035 }
00036     
00037 PCBlockDiagonal::~PCBlockDiagonal()
00038 {
00039     MatDestroy(mPCContext.A11_matrix_subblock);
00040     MatDestroy(mPCContext.A22_matrix_subblock);
00041     
00042     PCDestroy(mPCContext.PC_amg_A11);
00043     PCDestroy(mPCContext.PC_amg_A22);
00044 }
00045 
00046 void PCBlockDiagonal::PCBlockDiagonalCreate(KSP& rKspObject)
00047 {
00048     KSPGetPC(rKspObject, &mPetscPCObject);        
00049     
00050     Mat system_matrix, dummy;
00051     MatStructure flag;    
00052     KSPGetOperators(rKspObject, &system_matrix, &dummy, &flag); 
00053     
00054     PetscInt num_rows, num_columns;
00055     MatGetSize(system_matrix, &num_rows, &num_columns);
00056     assert(num_rows==num_columns);  
00057     
00058     PCSetType(mPetscPCObject, PCSHELL);
00059     
00060     // Register PC context so it gets passed to PCBlockDiagonalApply
00061     PCShellSetContext(mPetscPCObject, &mPCContext);
00062         
00063     // Get matrix sublock A11
00064     IS A11_rows, A11_columns;    
00065     ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 2, &A11_rows);
00066     ISCreateStride(PETSC_COMM_WORLD, num_columns/2, 0, 2, &A11_columns);
00067     
00068     MatGetSubMatrix(system_matrix, A11_rows, A11_columns, PETSC_DECIDE, MAT_INITIAL_MATRIX, &mPCContext.A11_matrix_subblock);
00069 
00070     ISDestroy(A11_rows);
00071     ISDestroy(A11_columns);    
00072 
00073     // Get matrix sublock A22
00074     IS A22_rows, A22_columns;    
00075     ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 1, 2, &A22_rows);
00076     ISCreateStride(PETSC_COMM_WORLD, num_columns/2, 1, 2, &A22_columns);
00077     
00078     MatGetSubMatrix(system_matrix, A22_rows, A22_columns, PETSC_DECIDE, MAT_INITIAL_MATRIX, &mPCContext.A22_matrix_subblock);
00079 
00080     ISDestroy(A22_rows);
00081     ISDestroy(A22_columns);    
00082 
00083     // Register call-back function
00084     PCShellSetApply(mPetscPCObject, PCBlockDiagonalApply);
00085 }
00086 
00087 void PCBlockDiagonal::PCBlockDiagonalSetUp()
00088 {
00089     // These options will get read by PCSetFromOptions
00090     PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");        
00091     PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");            
00092     PetscOptionsSetValue("-pc_hypre_type", "boomeramg");        
00093 
00094     // Set up amg preconditioner for block A11
00095     PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A11));
00096     PCSetType(mPCContext.PC_amg_A11, PCHYPRE);
00097     PCSetOperators(mPCContext.PC_amg_A11, mPCContext.A11_matrix_subblock, mPCContext.A11_matrix_subblock, DIFFERENT_NONZERO_PATTERN);//   SAME_PRECONDITIONER);
00098     PCSetFromOptions(mPCContext.PC_amg_A11);
00099     PCSetUp(mPCContext.PC_amg_A11);
00100 
00101     // Set up amg preconditioner for block A22
00102     PCCreate(PETSC_COMM_WORLD, &(mPCContext.PC_amg_A22));
00103     PCSetType(mPCContext.PC_amg_A22, PCHYPRE);
00104     PCSetOperators(mPCContext.PC_amg_A22, mPCContext.A22_matrix_subblock, mPCContext.A22_matrix_subblock, DIFFERENT_NONZERO_PATTERN);//   SAME_PRECONDITIONER);
00105     PCSetFromOptions(mPCContext.PC_amg_A22);
00106     PCSetUp(mPCContext.PC_amg_A22);        
00107 }
00108 
00109 PetscErrorCode PCBlockDiagonalApply(void *pc_context, Vec x, Vec y)
00110 {
00113     
00114     // Cast the pointer to a PC context to our defined type
00115     PCBlockDiagonal::PCBlockDiagonalContext* block_diag_context = (PCBlockDiagonal::PCBlockDiagonalContext*) pc_context;
00116     assert(block_diag_context!=NULL); 
00117     
00119     PetscInt num_rows;
00120     VecGetSize(x, &num_rows);
00121 
00122     Vec x11 = PetscTools::CreateVec(num_rows/2);
00123     Vec x22 = PetscTools::CreateVec(num_rows/2);
00124     Vec y11 = PetscTools::CreateVec(num_rows/2);
00125     Vec y22 = PetscTools::CreateVec(num_rows/2);
00126 
00127     IS A11_rows;    
00128     ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 2, &A11_rows);    
00129     
00130     VecScatter A11_scatter_ctx;
00131     VecScatterCreate(x, A11_rows, x11, PETSC_NULL, &A11_scatter_ctx);
00132 
00133 #if (PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)
00134     VecScatterBegin(A11_scatter_ctx, x, x11, INSERT_VALUES, SCATTER_FORWARD);
00135     VecScatterEnd(A11_scatter_ctx, x, x11, INSERT_VALUES, SCATTER_FORWARD);
00136 #else
00137     VecScatterBegin(x, x11, INSERT_VALUES, SCATTER_FORWARD, A11_scatter_ctx);
00138     VecScatterEnd(x, x11, INSERT_VALUES, SCATTER_FORWARD, A11_scatter_ctx);
00139 #endif    
00140 
00141     IS A22_rows;
00142     ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 1, 2, &A22_rows);
00143 
00144     VecScatter A22_scatter_ctx;
00145     VecScatterCreate(x, A22_rows, x22, PETSC_NULL, &A22_scatter_ctx);
00146 
00147 #if (PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)
00148     VecScatterBegin(A22_scatter_ctx, x, x22, INSERT_VALUES, SCATTER_FORWARD);
00149     VecScatterEnd(A22_scatter_ctx, x, x22, INSERT_VALUES, SCATTER_FORWARD);
00150 #else
00151     VecScatterBegin(x, x22, INSERT_VALUES, SCATTER_FORWARD, A22_scatter_ctx);
00152     VecScatterEnd(x, x22, INSERT_VALUES, SCATTER_FORWARD, A22_scatter_ctx);
00153 #endif    
00154     
00155     
00156     PCApply(block_diag_context->PC_amg_A11, x11, y11);
00157     PCApply(block_diag_context->PC_amg_A22, x22, y22);
00158 
00160 
00161 #if (PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)
00162     VecScatterBegin(A11_scatter_ctx, y11, y, INSERT_VALUES, SCATTER_REVERSE);
00163     VecScatterEnd(A11_scatter_ctx, y11, y, INSERT_VALUES, SCATTER_REVERSE);
00164 #else
00165     VecScatterBegin(y11, y, INSERT_VALUES, SCATTER_REVERSE, A11_scatter_ctx);
00166     VecScatterEnd(y11, y, INSERT_VALUES, SCATTER_REVERSE, A11_scatter_ctx);
00167 #endif    
00168 
00169 #if (PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)
00170     VecScatterBegin(A22_scatter_ctx, y22, y, INSERT_VALUES, SCATTER_REVERSE);
00171     VecScatterEnd(A22_scatter_ctx, y22, y, INSERT_VALUES, SCATTER_REVERSE);
00172 #else
00173     VecScatterBegin(y22, y, INSERT_VALUES, SCATTER_REVERSE, A22_scatter_ctx);
00174     VecScatterEnd(y22, y, INSERT_VALUES, SCATTER_REVERSE, A22_scatter_ctx);
00175 #endif    
00176     
00178 
00179     ISDestroy(A11_rows);
00180     ISDestroy(A22_rows);
00181         
00182     VecScatterDestroy(A11_scatter_ctx);
00183     VecScatterDestroy(A22_scatter_ctx);    
00184     
00185     VecDestroy(x11);
00186     VecDestroy(y11);
00187 
00188     VecDestroy(x22);
00189     VecDestroy(y22);
00190 
00191     
00192     return 0;
00193 }    

Generated on Tue Aug 4 16:10:22 2009 for Chaste by  doxygen 1.5.5