Chaste Release::3.1
|
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 STOKESFLOWSOLVER_HPP_ 00037 #define STOKESFLOWSOLVER_HPP_ 00038 00039 #include "AbstractContinuumMechanicsSolver.hpp" 00040 #include "LinearSystem.hpp" 00041 #include "LinearBasisFunction.hpp" 00042 #include "QuadraticBasisFunction.hpp" 00043 #include "Timer.hpp" 00044 #include "PetscMatTools.hpp" 00045 #include "PetscVecTools.hpp" 00046 #include "StokesFlowProblemDefinition.hpp" 00047 #include "StokesFlowAssembler.hpp" 00048 #include "StokesFlowPreconditionerAssembler.hpp" 00049 #include "ContinuumMechanicsNeumannBcsAssembler.hpp" 00050 00054 template<unsigned DIM> 00055 class StokesFlowSolver : public AbstractContinuumMechanicsSolver<DIM> 00056 { 00057 friend class TestStokesFlowSolver; 00058 00059 private: 00061 StokesFlowProblemDefinition<DIM>& mrProblemDefinition; 00062 00064 StokesFlowAssembler<DIM>* mpStokesFlowAssembler; 00065 00070 StokesFlowPreconditionerAssembler<DIM>* mpStokesFlowPreconditionerAssembler; 00071 00076 ContinuumMechanicsNeumannBcsAssembler<DIM>* mpNeumannBcsAssembler; 00077 00083 double mKspAbsoluteTol; 00084 00088 void AssembleSystem(); 00089 00090 public: 00091 00099 StokesFlowSolver(QuadraticMesh<DIM>& rQuadMesh, 00100 StokesFlowProblemDefinition<DIM>& rProblemDefinition, 00101 std::string outputDirectory); 00102 00106 virtual ~StokesFlowSolver(); 00107 00111 void Solve(); 00112 00119 void SetKspAbsoluteTolerance(double kspAbsoluteTolerance); 00120 00121 00126 std::vector<c_vector<double,DIM> >& rGetSpatialSolution(); 00127 00132 std::vector<c_vector<double,DIM> >& rGetVelocities(); 00133 }; 00134 00136 // Implementation 00138 00139 00140 template<unsigned DIM> 00141 StokesFlowSolver<DIM>::StokesFlowSolver(QuadraticMesh<DIM>& rQuadMesh, 00142 StokesFlowProblemDefinition<DIM>& rProblemDefinition, 00143 std::string outputDirectory) 00144 : AbstractContinuumMechanicsSolver<DIM>(rQuadMesh, rProblemDefinition, outputDirectory, INCOMPRESSIBLE), 00145 mrProblemDefinition(rProblemDefinition), 00146 mKspAbsoluteTol(-1) 00147 { 00148 assert(DIM==2 || DIM==3); 00149 assert(!mrProblemDefinition.rGetDirichletNodes().empty()); 00150 00151 mpStokesFlowAssembler = new StokesFlowAssembler<DIM>(&this->mrQuadMesh, &mrProblemDefinition); 00152 mpStokesFlowPreconditionerAssembler = new StokesFlowPreconditionerAssembler<DIM>(&this->mrQuadMesh, &mrProblemDefinition); 00153 mpNeumannBcsAssembler = new ContinuumMechanicsNeumannBcsAssembler<DIM>(&this->mrQuadMesh, &mrProblemDefinition); 00154 } 00155 00156 template<unsigned DIM> 00157 StokesFlowSolver<DIM>::~StokesFlowSolver() 00158 { 00159 delete mpStokesFlowAssembler; 00160 delete mpStokesFlowPreconditionerAssembler; 00161 delete mpNeumannBcsAssembler; 00162 } 00163 00164 template<unsigned DIM> 00165 void StokesFlowSolver<DIM>::Solve() 00166 { 00167 mrProblemDefinition.Validate(); 00168 00169 if(this->mVerbose) 00170 { 00171 Timer::Reset(); 00172 } 00173 00174 // Assemble Jacobian (and preconditioner) 00175 MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ASSEMBLE); 00176 AssembleSystem(); 00177 MechanicsEventHandler::EndEvent(MechanicsEventHandler::ASSEMBLE); 00178 00179 if(this->mVerbose) 00180 { 00181 Timer::PrintAndReset("AssembleSystem"); 00182 } 00183 00184 /* 00185 * Solve the linear system using PETSc GMRES and an LU factorisation 00186 * of the preconditioner. Note we don't call Solve on the linear_system 00187 * as we want to set PETSc options. 00188 */ 00189 MechanicsEventHandler::BeginEvent(MechanicsEventHandler::SOLVE); 00190 00191 Vec solution; 00192 VecDuplicate(this->mLinearSystemRhsVector,&solution); 00193 PetscVecTools::Zero(solution); 00194 00195 KSP solver; 00196 KSPCreate(PETSC_COMM_WORLD,&solver); 00197 KSPSetOperators(solver, this->mSystemLhsMatrix, this->mPreconditionMatrix, DIFFERENT_NONZERO_PATTERN /*in precond between successive solves*/); 00198 KSPSetType(solver, KSPGMRES); 00199 00200 if (mKspAbsoluteTol < 0) 00201 { 00202 double ksp_rel_tol = 1e-6; 00203 KSPSetTolerances(solver, ksp_rel_tol, PETSC_DEFAULT, PETSC_DEFAULT, 10000 /*max iter*/); //hopefully with the preconditioner this max is way too high 00204 } 00205 else 00206 { 00207 KSPSetTolerances(solver, 1e-16, mKspAbsoluteTol, PETSC_DEFAULT, 10000 /*max iter*/); //hopefully with the preconditioner this max is way too high 00208 } 00209 unsigned num_restarts = 100; 00210 KSPGMRESSetRestart(solver,num_restarts); // gmres num restarts 00211 00212 PC pc; 00213 KSPGetPC(solver, &pc); 00214 PCSetType(pc, PCJACOBI); 00215 00216 KSPSetUp(solver); 00217 00218 KSPSetFromOptions(solver); 00219 00220 // ///// For printing matrix when debugging 00221 // OutputFileHandler handler("TEMP"); 00222 // out_stream p_file = handler.OpenOutputFile("matrix.txt"); 00223 // for(unsigned i=0; i<this->mNumDofs; i++) 00224 // { 00225 // for(unsigned j=0; j<this->mNumDofs; j++) 00226 // { 00227 // *p_file << PetscMatTools::GetElement(this->mSystemLhsMatrix, i, j) << " "; 00228 // } 00229 // *p_file << "\n"; 00230 // } 00231 // p_file->close(); 00232 // 00233 // out_stream p_file2 = handler.OpenOutputFile("rhs.txt"); 00234 // for(unsigned i=0; i<this->mNumDofs; i++) 00235 // { 00236 // *p_file2 << PetscVecTools::GetElement(this->mLinearSystemRhsVector, i) << "\n"; 00237 // } 00238 // p_file2->close(); 00239 00240 if(this->mVerbose) 00241 { 00242 Timer::PrintAndReset("KSP Setup"); 00243 } 00244 00245 KSPSolve(solver,this->mLinearSystemRhsVector,solution); 00246 00247 KSPConvergedReason reason; 00248 KSPGetConvergedReason(solver,&reason); 00249 KSPEXCEPT(reason); 00250 00251 if(this->mVerbose) 00252 { 00253 Timer::PrintAndReset("KSP Solve"); 00254 int num_iters; 00255 KSPGetIterationNumber(solver, &num_iters); 00256 std::cout << "[" << PetscTools::GetMyRank() << "]: Num iterations = " << num_iters << "\n" << std::flush; 00257 } 00258 00259 MechanicsEventHandler::EndEvent(MechanicsEventHandler::SOLVE); 00260 00261 // Copy solution into the std::vector 00262 ReplicatableVector solution_repl(solution); 00263 for (unsigned i=0; i<this->mNumDofs; i++) 00264 { 00265 this->mCurrentSolution[i] = solution_repl[i]; 00266 } 00267 00268 // Remove pressure dummy values (P=0 at internal nodes, which should have been 00269 // been the result of the solve above), by linear interpolating from vertices of 00270 // edges to the internal node 00271 this->RemovePressureDummyValuesThroughLinearInterpolation(); 00272 00273 PetscTools::Destroy(solution); 00274 KSPDestroy(PETSC_DESTROY_PARAM(solver)); 00275 00276 this->WriteCurrentSpatialSolution("flow_solution", "nodes"); 00277 this->WriteCurrentPressureSolution(); 00278 } 00279 00280 00281 00282 template<unsigned DIM> 00283 void StokesFlowSolver<DIM>::AssembleSystem() 00284 { 00285 // Use assembler to assemble volume integral part.... 00286 mpStokesFlowAssembler->SetMatrixToAssemble(this->mSystemLhsMatrix, true); 00287 mpStokesFlowAssembler->SetVectorToAssemble(this->mLinearSystemRhsVector, true); 00288 mpStokesFlowAssembler->Assemble(); 00289 00290 mpStokesFlowPreconditionerAssembler->SetMatrixToAssemble(this->mPreconditionMatrix, true); 00291 mpStokesFlowPreconditionerAssembler->AssembleMatrix(); 00292 00293 mpNeumannBcsAssembler->SetVectorToAssemble(this->mLinearSystemRhsVector, false /*don't zero!*/); 00294 mpNeumannBcsAssembler->AssembleVector(); 00295 00296 PetscVecTools::Finalise(this->mLinearSystemRhsVector); 00297 PetscMatTools::SwitchWriteMode(this->mSystemLhsMatrix); 00298 PetscMatTools::SwitchWriteMode(this->mPreconditionMatrix); 00299 00300 // Note: maintaining symmetry for Dirichlet BCs is possible at the moment (the second parameter) 00301 // but doing so for the identity block is not yet implemented (the second parameter must be false) 00302 // Not sure if maintaining symmetry is worth it - may allow CG to work, but matrix is indefinite 00303 // so GC not guaranteed to work.. 00304 // 00305 // Note: the identity block needs to be added before the BCs - see comments in 00306 // PetscMatTools::ZeroRowsWithValueOnDiagonal() 00307 this->AddIdentityBlockForDummyPressureVariables(LINEAR_PROBLEM); 00308 this->ApplyDirichletBoundaryConditions(LINEAR_PROBLEM, true); 00309 00310 PetscVecTools::Finalise(this->mLinearSystemRhsVector); 00311 PetscMatTools::Finalise(this->mSystemLhsMatrix); 00312 PetscMatTools::Finalise(this->mPreconditionMatrix); 00313 } 00314 00315 template<unsigned DIM> 00316 void StokesFlowSolver<DIM>::SetKspAbsoluteTolerance(double kspAbsoluteTolerance) 00317 { 00318 assert(kspAbsoluteTolerance > 0); 00319 mKspAbsoluteTol = kspAbsoluteTolerance; 00320 } 00321 00322 template<unsigned DIM> 00323 std::vector<c_vector<double,DIM> >& StokesFlowSolver<DIM>::rGetSpatialSolution() 00324 { 00325 this->mSpatialSolution.resize(this->mrQuadMesh.GetNumNodes(), zero_vector<double>(DIM)); 00326 for (unsigned i=0; i<this->mrQuadMesh.GetNumNodes(); i++) 00327 { 00328 for (unsigned j=0; j<DIM; j++) 00329 { 00330 // DIM+1 is the problem dimension 00331 this->mSpatialSolution[i](j) = this->mCurrentSolution[(DIM+1)*i+j]; 00332 } 00333 } 00334 return this->mSpatialSolution; 00335 } 00336 00337 template<unsigned DIM> 00338 std::vector<c_vector<double,DIM> >& StokesFlowSolver<DIM>::rGetVelocities() 00339 { 00340 return rGetSpatialSolution(); 00341 } 00342 00343 00344 00345 #endif /* STOKESFLOWSOLVER_HPP_ */