Chaste Release::3.1
AbstractContinuumMechanicsSolver.hpp
00001 /*
00002 
00003 Copyright (c) 2005-2012, University of Oxford.
00004 All rights reserved.
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Redistribution and use in source and binary forms, with or without
00013 modification, are permitted provided that the following conditions are met:
00014  * Redistributions of source code must retain the above copyright notice,
00015    this list of conditions and the following disclaimer.
00016  * Redistributions in binary form must reproduce the above copyright notice,
00017    this list of conditions and the following disclaimer in the documentation
00018    and/or other materials provided with the distribution.
00019  * Neither the name of the University of Oxford nor the names of its
00020    contributors may be used to endorse or promote products derived from this
00021    software without specific prior written permission.
00022 
00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00033 
00034 */
00035 
00036 #ifndef ABSTRACTCONTINUUMMECHANICSSOLVER_HPP_
00037 #define ABSTRACTCONTINUUMMECHANICSSOLVER_HPP_
00038 
00039 #include "ContinuumMechanicsProblemDefinition.hpp"
00040 #include "CompressibilityType.hpp"
00041 #include "LinearSystem.hpp"
00042 #include "OutputFileHandler.hpp"
00043 #include "ReplicatableVector.hpp"
00044 #include "QuadraticMesh.hpp"
00045 #include "Warnings.hpp"
00046 #include "PetscException.hpp"
00047 #include "GaussianQuadratureRule.hpp"
00048 #include "PetscTools.hpp"
00049 #include "MechanicsEventHandler.hpp"
00050 #include "CommandLineArguments.hpp"
00051 #include "VtkMeshWriter.hpp"
00052 
00053 
00059 typedef enum _ApplyDirichletBcsType
00060 {
00061     LINEAR_PROBLEM,
00062     NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY,
00063     NONLINEAR_PROBLEM_APPLY_TO_EVERYTHING,
00064 } ApplyDirichletBcsType;
00065 
00066 
00081 template<unsigned DIM>
00082 class AbstractContinuumMechanicsSolver
00083 {
00084 protected:
00089     QuadraticMesh<DIM>& mrQuadMesh;
00090 
00092     ContinuumMechanicsProblemDefinition<DIM>& mrProblemDefinition;
00093 
00095     bool mWriteOutput;
00096 
00098     std::string mOutputDirectory;
00099 
00101     OutputFileHandler* mpOutputFileHandler;
00102 
00106     std::vector<c_vector<double,DIM> > mSpatialSolution;
00107 
00109     std::vector<double> mPressureSolution;
00110 
00111 
00118     std::vector<double> mCurrentSolution;
00119 
00121     GaussianQuadratureRule<DIM>* mpQuadratureRule;
00122 
00124     GaussianQuadratureRule<DIM-1>* mpBoundaryQuadratureRule;
00125 
00130     CompressibilityType mCompressibilityType;
00131 
00140     unsigned mProblemDimension;
00141 
00145     unsigned mNumDofs;
00146 
00152     bool mVerbose;
00153 
00174     Vec mResidualVector;
00175 
00179     Vec mLinearSystemRhsVector;
00180 
00184     Mat mSystemLhsMatrix;
00185 
00189     Vec mDirichletBoundaryConditionsVector;
00190 
00194     Mat mPreconditionMatrix;
00195 
00199     void AllocateMatrixMemory();
00200 
00201 
00227     void ApplyDirichletBoundaryConditions(ApplyDirichletBcsType type, bool symmetricProblem);
00228 
00261     void AddIdentityBlockForDummyPressureVariables(ApplyDirichletBcsType type);
00262 
00263 
00280     void RemovePressureDummyValuesThroughLinearInterpolation();
00281 
00282 public:
00290     AbstractContinuumMechanicsSolver(QuadraticMesh<DIM>& rQuadMesh,
00291                                      ContinuumMechanicsProblemDefinition<DIM>& rProblemDefinition,
00292                                      std::string outputDirectory,
00293                                      CompressibilityType compressibilityType);
00294 
00296     virtual ~AbstractContinuumMechanicsSolver();
00297 
00309     void WriteCurrentSpatialSolution(std::string fileName, std::string fileExtension, int counterToAppend=-1);
00310 
00323     void WriteCurrentPressureSolution( int counterToAppend=-1);
00324 
00330     void SetWriteOutput(bool writeOutput=true);
00331 
00336     void CreateVtkOutput(std::string spatialSolutionName="Spatial solution");
00337 
00342     std::vector<double>& rGetCurrentSolution()
00343     {
00344         return mCurrentSolution;
00345     }
00346 
00351     virtual std::vector<c_vector<double,DIM> >& rGetSpatialSolution()=0;
00352 
00361     std::vector<double>& rGetPressures();
00362 
00363 };
00364 
00365 
00366 template<unsigned DIM>
00367 AbstractContinuumMechanicsSolver<DIM>::AbstractContinuumMechanicsSolver(QuadraticMesh<DIM>& rQuadMesh,
00368                                                                         ContinuumMechanicsProblemDefinition<DIM>& rProblemDefinition,
00369                                                                         std::string outputDirectory,
00370                                                                         CompressibilityType compressibilityType)
00371     : mrQuadMesh(rQuadMesh),
00372       mrProblemDefinition(rProblemDefinition),
00373       mOutputDirectory(outputDirectory),
00374       mpOutputFileHandler(NULL),
00375       mpQuadratureRule(NULL),
00376       mpBoundaryQuadratureRule(NULL),
00377       mCompressibilityType(compressibilityType),
00378       mResidualVector(NULL),
00379       mSystemLhsMatrix(NULL),
00380       mPreconditionMatrix(NULL)
00381 {
00382     assert(DIM==2 || DIM==3);
00383 
00384     mVerbose = (mrProblemDefinition.GetVerboseDuringSolve() ||
00385                 CommandLineArguments::Instance()->OptionExists("-mech_verbose") ||
00386                 CommandLineArguments::Instance()->OptionExists("-mech_very_verbose") );
00387 
00388     mWriteOutput = (mOutputDirectory != "");
00389     if (mWriteOutput)
00390     {
00391         mpOutputFileHandler = new OutputFileHandler(mOutputDirectory);
00392     }
00393 
00394     // see dox for mProblemDimension
00395     mProblemDimension = mCompressibilityType==COMPRESSIBLE ? DIM : DIM+1;
00396     mNumDofs = mProblemDimension*mrQuadMesh.GetNumNodes();
00397 
00398     AllocateMatrixMemory();
00399 
00400     mpQuadratureRule = new GaussianQuadratureRule<DIM>(3);
00401     mpBoundaryQuadratureRule = new GaussianQuadratureRule<DIM-1>(3);
00402 
00403     mCurrentSolution.resize(mNumDofs, 0.0);
00404 }
00405 
00406 
00407 template<unsigned DIM>
00408 AbstractContinuumMechanicsSolver<DIM>::~AbstractContinuumMechanicsSolver()
00409 {
00410     if (mpOutputFileHandler)
00411     {
00412         delete mpOutputFileHandler;
00413     }
00414 
00415     if (mpQuadratureRule)
00416     {
00417         delete mpQuadratureRule;
00418         delete mpBoundaryQuadratureRule;
00419     }
00420 
00421     if (mResidualVector)
00422     {
00423         PetscTools::Destroy(mResidualVector);
00424         PetscTools::Destroy(mLinearSystemRhsVector);
00425         PetscTools::Destroy(mSystemLhsMatrix);
00426         PetscTools::Destroy(mPreconditionMatrix);
00427     }
00428 
00429     if (mDirichletBoundaryConditionsVector)
00430     {
00431         PetscTools::Destroy(mDirichletBoundaryConditionsVector);
00432     }
00433 }
00434 
00435 template<unsigned DIM>
00436 void AbstractContinuumMechanicsSolver<DIM>::WriteCurrentSpatialSolution(std::string fileName,
00437                                                                         std::string fileExtension,
00438                                                                         int counterToAppend)
00439 {
00440     // Only write output if the flag mWriteOutput has been set
00441     if (!mWriteOutput)
00442     {
00443         return;
00444     }
00445 
00446     if (PetscTools::AmMaster())
00447     {
00448         std::stringstream file_name;
00449 
00450         file_name << fileName;
00451         if (counterToAppend >= 0)
00452         {
00453             file_name << "_" << counterToAppend;
00454         }
00455         file_name << "." << fileExtension;
00456 
00457         out_stream p_file = mpOutputFileHandler->OpenOutputFile(file_name.str());
00458 
00459         std::vector<c_vector<double,DIM> >& r_spatial_solution = rGetSpatialSolution();
00460         for (unsigned i=0; i<r_spatial_solution.size(); i++)
00461         {
00462     //        for (unsigned j=0; j<DIM; j++)
00463     //        {
00464     //            *p_file << mrQuadMesh.GetNode(i)->rGetLocation()[j] << " ";
00465     //        }
00466             for (unsigned j=0; j<DIM; j++)
00467             {
00468                 *p_file << r_spatial_solution[i](j) << " ";
00469             }
00470             *p_file << "\n";
00471         }
00472         p_file->close();
00473     }
00474     PetscTools::Barrier("WriteSpatial");
00475 }
00476 
00477 template<unsigned DIM>
00478 void AbstractContinuumMechanicsSolver<DIM>::WriteCurrentPressureSolution(int counterToAppend)
00479 {
00480     // Only write output if the flag mWriteOutput has been set
00481     if (!mWriteOutput)
00482     {
00483         return;
00484     }
00485 
00486     if (PetscTools::AmMaster())
00487     {
00488         std::stringstream file_name;
00489 
00490         file_name << "pressure";
00491         if (counterToAppend >= 0)
00492         {
00493             file_name << "_" << counterToAppend;
00494         }
00495         file_name << ".txt";
00496 
00497         out_stream p_file = mpOutputFileHandler->OpenOutputFile(file_name.str());
00498 
00499         std::vector<double>& r_pressure = rGetPressures();
00500         for (unsigned i=0; i<r_pressure.size(); i++)
00501         {
00502             for (unsigned j=0; j<DIM; j++)
00503             {
00504                 *p_file << mrQuadMesh.GetNode(i)->rGetLocation()[j] << " ";
00505             }
00506 
00507             *p_file << r_pressure[i] << "\n";
00508         }
00509         p_file->close();
00510     }
00511     PetscTools::Barrier("WritePressure");
00512 }
00513 template<unsigned DIM>
00514 void AbstractContinuumMechanicsSolver<DIM>::SetWriteOutput(bool writeOutput)
00515 {
00516     if (writeOutput && (mOutputDirectory==""))
00517     {
00518         EXCEPTION("Can't write output if no output directory was given in constructor");
00519     }
00520     mWriteOutput = writeOutput;
00521 }
00522 
00523 template<unsigned DIM>
00524 void AbstractContinuumMechanicsSolver<DIM>::CreateVtkOutput(std::string spatialSolutionName)
00525 {
00526     if (this->mOutputDirectory=="")
00527     {
00528         EXCEPTION("No output directory was given so no output was written, cannot convert to VTK format");
00529     }
00530 #ifdef CHASTE_VTK
00531     VtkMeshWriter<DIM, DIM> mesh_writer(this->mOutputDirectory + "/vtk", "solution", true);
00532 
00533     mesh_writer.AddPointData(spatialSolutionName, this->rGetSpatialSolution());
00534 
00535     if (mCompressibilityType==INCOMPRESSIBLE)
00536     {
00537         mesh_writer.AddPointData("Pressure", rGetPressures());
00538     }
00539 
00540     //Output the element attribute as cell data.
00541     std::vector<double> element_attribute;
00542     for(typename QuadraticMesh<DIM>::ElementIterator iter = this->mrQuadMesh.GetElementIteratorBegin();
00543         iter != this->mrQuadMesh.GetElementIteratorEnd();
00544         ++iter)
00545     {
00546         element_attribute.push_back(iter->GetAttribute());
00547     }
00548     mesh_writer.AddCellData("Attribute", element_attribute);
00549 
00550     mesh_writer.WriteFilesUsingMesh(this->mrQuadMesh);
00551 #endif
00552 }
00553 
00554 
00555 template<unsigned DIM>
00556 std::vector<double>& AbstractContinuumMechanicsSolver<DIM>::rGetPressures()
00557 {
00558     assert(mProblemDimension==DIM+1);
00559 
00560     mPressureSolution.clear();
00561     mPressureSolution.resize(mrQuadMesh.GetNumNodes());
00562 
00563     for (unsigned i=0; i<mrQuadMesh.GetNumNodes(); i++)
00564     {
00565         mPressureSolution[i] = mCurrentSolution[mProblemDimension*i + DIM];
00566     }
00567     return mPressureSolution;
00568 }
00569 
00570 
00571 template<unsigned DIM>
00572 void AbstractContinuumMechanicsSolver<DIM>::RemovePressureDummyValuesThroughLinearInterpolation()
00573 {
00574     assert(mProblemDimension==DIM+1);
00575 
00576     // For quadratic triangles, node 3 is between nodes 1 and 2, node 4 is between 0 and 2, etc
00577     unsigned internal_nodes_2d[3] = {3,4,5};
00578     unsigned neighbouring_vertices_2d[3][2] = { {1,2}, {2,0}, {0,1} };
00579 
00580     // ordering for quadratic tetrahedra
00581     unsigned internal_nodes_3d[6] = {4,5,6,7,8,9};
00582     unsigned neighbouring_vertices_3d[6][2] = { {0,1}, {1,2}, {0,2}, {0,3}, {1,3}, {2,3} };
00583 
00584     unsigned num_internal_nodes_per_element = DIM==2 ? 3 : 6;
00585 
00586     // loop over elements, then loop over edges.
00587     for (typename AbstractTetrahedralMesh<DIM,DIM>::ElementIterator iter = mrQuadMesh.GetElementIteratorBegin();
00588          iter != mrQuadMesh.GetElementIteratorEnd();
00589          ++iter)
00590     {
00591         for(unsigned i=0; i<num_internal_nodes_per_element; i++)
00592         {
00593             unsigned global_index;
00594             double left_val;
00595             double right_val;
00596 
00597             if(DIM==2)
00598             {
00599                 global_index = iter->GetNodeGlobalIndex( internal_nodes_2d[i] );
00600                 unsigned vertex_0_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_2d[i][0] );
00601                 unsigned vertex_1_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_2d[i][1] );
00602                 left_val = mCurrentSolution[mProblemDimension*vertex_0_global_index + DIM];
00603                 right_val = mCurrentSolution[mProblemDimension*vertex_1_global_index + DIM];
00604             }
00605             else
00606             {
00607                 global_index = iter->GetNodeGlobalIndex( internal_nodes_3d[i] );
00608                 unsigned vertex_0_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_3d[i][0] );
00609                 unsigned vertex_1_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_3d[i][1] );
00610                 left_val = mCurrentSolution[mProblemDimension*vertex_0_global_index + DIM];
00611                 right_val = mCurrentSolution[mProblemDimension*vertex_1_global_index + DIM];
00612             }
00613 
00614             // this line assumes the internal node is midway between the two vertices
00615             mCurrentSolution[mProblemDimension*global_index + DIM] =  0.5 * (left_val + right_val);
00616         }
00617     }
00618 }
00619 
00620 /*
00621  * This method applies the appropriate BCs to the residual and/or linear system
00622  *
00623  * For the latter, and second input parameter==true, the BCs are imposed in such a way as
00624  * to ensure that a symmetric linear system remains symmetric. For each row with boundary condition
00625  * applied, both the row and column are zero'd and the RHS vector modified to take into account the
00626  * zero'd column.
00627  *
00628  * Suppose we have a matrix
00629  * [a b c] [x] = [ b1 ]
00630  * [d e f] [y]   [ b2 ]
00631  * [g h i] [z]   [ b3 ]
00632  * and we want to apply the boundary condition x=v without losing symmetry if the matrix is
00633  * symmetric. We apply the boundary condition
00634  * [1 0 0] [x] = [ v  ]
00635  * [d e f] [y]   [ b2 ]
00636  * [g h i] [z]   [ b3 ]
00637  * and then zero the column as well, adding a term to the RHS to take account for the
00638  * zero-matrix components
00639  * [1 0 0] [x] = [ v  ] - v[ 0 ]
00640  * [0 e f] [y]   [ b2 ]    [ d ]
00641  * [0 h i] [z]   [ b3 ]    [ g ]
00642  * Note the last term is the first column of the matrix, with one component zeroed, and
00643  * multiplied by the boundary condition value. This last term is then stored in
00644  * mDirichletBoundaryConditionsVector, and in general is equal to:
00645  * SUM_{d=1..D} v_d a'_d
00646  * where v_d is the boundary value of boundary condition d (d an index into the matrix),
00647  * and a'_d is the dth-column of the matrix but with the d-th component zeroed, and where
00648  * there are D boundary conditions
00649  */
00650 template<unsigned DIM>
00651 void AbstractContinuumMechanicsSolver<DIM>::ApplyDirichletBoundaryConditions(ApplyDirichletBcsType type, bool symmetricProblem)
00652 {
00653     std::vector<unsigned> rows;
00654     std::vector<double> values;
00655 
00656     // Whether to apply symmetrically, ie alter columns as well as rows (see comment above)
00657     bool applySymmetrically = (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY) && symmetricProblem;
00658 
00659     if (applySymmetrically)
00660     {
00661         if(mDirichletBoundaryConditionsVector==NULL)
00662         {
00663             VecDuplicate(mResidualVector, &mDirichletBoundaryConditionsVector);
00664         }
00665 
00666         PetscVecTools::Zero(mDirichletBoundaryConditionsVector);
00667         PetscMatTools::Finalise(mSystemLhsMatrix);
00668     }
00669 
00671     // collect the entries to be altered
00673 
00674     for (unsigned i=0; i<mrProblemDefinition.rGetDirichletNodes().size(); i++)
00675     {
00676         unsigned node_index = mrProblemDefinition.rGetDirichletNodes()[i];
00677 
00678         for (unsigned j=0; j<DIM; j++)
00679         {
00680             double dirichlet_val = mrProblemDefinition.rGetDirichletNodeValues()[i](j);
00681 
00682             if(dirichlet_val != ContinuumMechanicsProblemDefinition<DIM>::FREE)
00683             {
00684                 double val;
00685                 unsigned dof_index = mProblemDimension*node_index+j;
00686 
00687                 if(type == LINEAR_PROBLEM)
00688                 {
00689                     val = dirichlet_val;
00690                 }
00691                 else
00692                 {
00693                     // The boundary conditions on the NONLINEAR SYSTEM are x=boundary_values
00694                     // on the boundary nodes. However:
00695                     // The boundary conditions on the LINEAR SYSTEM at each Newton step (Ju=f,
00696                     // where J is the Jacobian, u the negative update vector and f is the residual) is:
00697                     // u=current_soln-boundary_values on the boundary nodes
00698                     val = mCurrentSolution[dof_index] - dirichlet_val;
00699                 }
00700                 rows.push_back(dof_index);
00701                 values.push_back(val);
00702             }
00703         }
00704     }
00705 
00707     // do the alterations
00709 
00710     if (applySymmetrically)
00711     {
00712         // Modify the matrix columns
00713         for (unsigned i=0; i<rows.size(); i++)
00714         {
00715             unsigned col = rows[i];
00716             double minus_value = -values[i];
00717 
00718             // Get a vector which will store the column of the matrix (column d, where d is
00719             // the index of the row (and column) to be altered for the boundary condition.
00720             // Since the matrix is symmetric when get row number "col" and treat it as a column.
00721             // PETSc uses compressed row format and therefore getting rows is far more efficient
00722             // than getting columns.
00723             Vec matrix_col = PetscMatTools::GetMatrixRowDistributed(mSystemLhsMatrix,col);
00724 
00725             // Zero the correct entry of the column
00726             PetscVecTools::SetElement(matrix_col, col, 0.0);
00727 
00728             // Set up the RHS Dirichlet boundary conditions vector
00729             // Eg, for a boundary node at the zeroth node (x_0 = value), this is equal to
00730             //   -value*[0 a_21 a_31 .. a_N1]
00731             // and will be added to the RHS.
00732             PetscVecTools::AddScaledVector(mDirichletBoundaryConditionsVector, matrix_col, minus_value);
00733             PetscTools::Destroy(matrix_col);
00734         }
00735     }
00736 
00737     if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY) // ie doing a whole linear system
00738     {
00739         // Now zero the appropriate rows and columns of the matrix. If the matrix is symmetric we apply the
00740         // boundary conditions in a way the symmetry isn't lost (rows and columns). If not only the row is
00741         // zeroed.
00742         if (applySymmetrically)
00743         {
00744             PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal(mSystemLhsMatrix, rows, 1.0);
00745             PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal(mPreconditionMatrix, rows, 1.0);
00746 
00747             // Apply the RHS boundary conditions modification if required.
00748             PetscVecTools::AddScaledVector(mLinearSystemRhsVector, mDirichletBoundaryConditionsVector, 1.0);
00749         }
00750         else
00751         {
00752             PetscMatTools::ZeroRowsWithValueOnDiagonal(mSystemLhsMatrix, rows, 1.0);
00753             PetscMatTools::ZeroRowsWithValueOnDiagonal(mPreconditionMatrix, rows, 1.0);
00754         }
00755     }
00756 
00757     if (type!=LINEAR_PROBLEM)
00758     {
00759         for (unsigned i=0; i<rows.size(); i++)
00760         {
00761             PetscVecTools::SetElement(mResidualVector, rows[i], values[i]);
00762         }
00763     }
00764 
00765     if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
00766     {
00767         for (unsigned i=0; i<rows.size(); i++)
00768         {
00769             PetscVecTools::SetElement(mLinearSystemRhsVector, rows[i], values[i]);
00770         }
00771     }
00772 }
00773 
00774 
00775 template<unsigned DIM>
00776 void AbstractContinuumMechanicsSolver<DIM>::AddIdentityBlockForDummyPressureVariables(ApplyDirichletBcsType type)
00777 {
00778     assert(mCompressibilityType==INCOMPRESSIBLE);
00779 
00780     int lo, hi;
00781     VecGetOwnershipRange(mResidualVector, &lo, &hi);
00782 
00783     for(unsigned i=0; i<mrQuadMesh.GetNumNodes(); i++)
00784     {
00785         if(mrQuadMesh.GetNode(i)->IsInternal())
00786         {
00787             unsigned row = (DIM+1)*i + DIM; // DIM+1 is the problem dimension
00788             if(lo <= (int)row && (int)row < hi)
00789             {
00790                 if(type!=LINEAR_PROBLEM)
00791                 {
00792                     PetscVecTools::SetElement(mResidualVector, row, mCurrentSolution[row]-0.0);
00793                 }
00794                 if(type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY) // ie doing a whole linear system
00795                 {
00796                     double rhs_vector_val = type==LINEAR_PROBLEM ? 0.0 : mCurrentSolution[row]-0.0;
00797                     PetscVecTools::SetElement(mLinearSystemRhsVector, row, rhs_vector_val);
00798                     // this assumes the row is already zero, which is should be..
00799                     PetscMatTools::SetElement(mSystemLhsMatrix, row, row, 1.0);
00800                     PetscMatTools::SetElement(mPreconditionMatrix, row, row, 1.0);
00801                 }
00802             }
00803         }
00804     }
00805 }
00806 
00807 
00808 
00809 template<unsigned DIM>
00810 void AbstractContinuumMechanicsSolver<DIM>::AllocateMatrixMemory()
00811 {
00812     Vec template_vec = mrQuadMesh.GetDistributedVectorFactory()->CreateVec(mProblemDimension);
00813 
00815     // three vectors
00817     VecDuplicate(template_vec, &mResidualVector);
00818     VecDuplicate(mResidualVector, &mLinearSystemRhsVector);
00819     // the one is only allocated if it will be needed (in ApplyDirichletBoundaryConditions),
00820     // depending on whether the matrix is kept symmetric.
00821     mDirichletBoundaryConditionsVector = NULL;
00822     PetscTools::Destroy(template_vec);
00823 
00825     // two matrices
00827 
00828     int lo, hi;
00829     VecGetOwnershipRange(mResidualVector, &lo, &hi);
00830     PetscInt local_size = hi - lo;
00831 
00832 
00833     if (DIM==2)
00834     {
00835         // 2D: N elements around a point => 7N+3 non-zeros in that row? Assume N<=10 (structured mesh would have N_max=6) => 73.
00836         unsigned num_non_zeros = std::min(75u, mNumDofs);
00837 
00838         PetscTools::SetupMat(mSystemLhsMatrix, mNumDofs, mNumDofs, num_non_zeros, local_size, local_size);
00839         PetscTools::SetupMat(mPreconditionMatrix, mNumDofs, mNumDofs, num_non_zeros, local_size, local_size);
00840     }
00841     else
00842     {
00843         assert(DIM==3);
00844 
00845         // in 3d we get the number of containing elements for each node and use that to obtain an upper bound
00846         // for the number of non-zeros for each DOF associated with that node.
00847 
00848         int* num_non_zeros_each_row = new int[mNumDofs];
00849         for (unsigned i=0; i<mNumDofs; i++)
00850         {
00851             num_non_zeros_each_row[i] = 0;
00852         }
00853 
00854         for (unsigned i=0; i<mrQuadMesh.GetNumNodes(); i++)
00855         {
00856             // this upper bound neglects the fact that two containing elements will share the same nodes..
00857             // 4 = max num dofs associated with this node
00858             // 30 = 3*9+3 = 3 dimensions x 9 other nodes on this element   +  3 vertices with a pressure unknown
00859             unsigned num_non_zeros_upper_bound = 4 + 30*mrQuadMesh.GetNode(i)->GetNumContainingElements();
00860 
00861             num_non_zeros_upper_bound = std::min(num_non_zeros_upper_bound, mNumDofs);
00862 
00863             num_non_zeros_each_row[mProblemDimension*i + 0] = num_non_zeros_upper_bound;
00864             num_non_zeros_each_row[mProblemDimension*i + 1] = num_non_zeros_upper_bound;
00865             num_non_zeros_each_row[mProblemDimension*i + 2] = num_non_zeros_upper_bound;
00866 
00867             if (mCompressibilityType==INCOMPRESSIBLE)
00868             {
00869                 if(!mrQuadMesh.GetNode(i)->IsInternal())
00870                 {
00871                     num_non_zeros_each_row[mProblemDimension*i + 3] = num_non_zeros_upper_bound;
00872                 }
00873                 else
00874                 {
00875                     num_non_zeros_each_row[mProblemDimension*i + 3] = 1;
00876                 }
00877             }
00878         }
00879 
00880         // NOTE: PetscTools::SetupMat() or the below creates a MATAIJ matrix, which means the matrix will
00881         // be of type MATSEQAIJ if num_procs=1 and MATMPIAIJ otherwise. In the former case
00882         // MatSeqAIJSetPreallocation MUST be called [MatMPIAIJSetPreallocation will have
00883         // no effect (silently)], and vice versa in the latter case
00884 
00887         //PetscTools::SetupMat(mSystemLhsMatrix, mNumDofs, mNumDofs, 0, PETSC_DECIDE, PETSC_DECIDE);
00888         //PetscTools::SetupMat(mPreconditionMatrix, mNumDofs, mNumDofs, 0, PETSC_DECIDE, PETSC_DECIDE);
00890 
00891         // possible todo: create a PetscTools::SetupMatNoAllocation()
00892 
00893         // In the future, when parallelising, remember to think about MAT_IGNORE_OFF_PROC_ENTRIES (see #1682)
00894 
00895 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00896         MatCreate(PETSC_COMM_WORLD,local_size,local_size,mNumDofs,mNumDofs,&mSystemLhsMatrix);
00897         MatCreate(PETSC_COMM_WORLD,local_size,local_size,mNumDofs,mNumDofs,&mPreconditionMatrix);
00898 #else //New API
00899         MatCreate(PETSC_COMM_WORLD,&mSystemLhsMatrix);
00900         MatCreate(PETSC_COMM_WORLD,&mPreconditionMatrix);
00901         MatSetSizes(mSystemLhsMatrix,local_size,local_size,mNumDofs,mNumDofs);
00902         MatSetSizes(mPreconditionMatrix,local_size,local_size,mNumDofs,mNumDofs);
00903 #endif
00904 
00905         if (PetscTools::IsSequential())
00906         {
00907             MatSetType(mSystemLhsMatrix, MATSEQAIJ);
00908             MatSetType(mPreconditionMatrix, MATSEQAIJ);
00909             MatSeqAIJSetPreallocation(mSystemLhsMatrix,    PETSC_NULL, num_non_zeros_each_row);
00910             MatSeqAIJSetPreallocation(mPreconditionMatrix, PETSC_NULL, num_non_zeros_each_row);
00911         }
00912         else
00913         {
00914             int* num_non_zeros_each_row_this_proc = new int[local_size];
00915             int* zero = new int[local_size];
00916             for (unsigned i=0; i<unsigned(local_size); i++)
00917             {
00918                 num_non_zeros_each_row_this_proc[i] = num_non_zeros_each_row[lo+i];
00919                 if(num_non_zeros_each_row_this_proc[i] > local_size)
00920                 {
00921                     num_non_zeros_each_row_this_proc[i] = local_size;
00922                 }
00923                 zero[i] = 0;
00924             }
00925 
00926             MatSetType(mSystemLhsMatrix, MATMPIAIJ);
00927             MatSetType(mPreconditionMatrix, MATMPIAIJ);
00928             MatMPIAIJSetPreallocation(mSystemLhsMatrix,    PETSC_NULL, num_non_zeros_each_row_this_proc, PETSC_NULL, num_non_zeros_each_row_this_proc);
00929             MatMPIAIJSetPreallocation(mPreconditionMatrix, PETSC_NULL, num_non_zeros_each_row_this_proc, PETSC_NULL, num_non_zeros_each_row_this_proc);
00930         }
00931 
00932         MatSetFromOptions(mSystemLhsMatrix);
00933         MatSetFromOptions(mPreconditionMatrix);
00934 
00935         //unsigned total_non_zeros = 0;
00936         //for (unsigned i=0; i<mNumDofs; i++)
00937         //{
00938         //   total_non_zeros += num_non_zeros_each_row[i];
00939         //}
00940         //std::cout << total_non_zeros << " versus " << 500*mNumDofs << "\n" << std::flush;
00941 
00942         delete [] num_non_zeros_each_row;
00943     }
00944 }
00945 #endif // ABSTRACTCONTINUUMMECHANICSSOLVER_HPP_