Chaste Release::3.1
BoundaryConditionsContainerImplementation.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 _BOUNDARYCONDITIONSCONTAINERIMPLEMENTATION_HPP_
00037 #define _BOUNDARYCONDITIONSCONTAINERIMPLEMENTATION_HPP_
00038 
00039 #include "PetscVecTools.hpp"
00040 #include "PetscMatTools.hpp"
00041 #include "BoundaryConditionsContainer.hpp"
00042 #include "DistributedVector.hpp"
00043 #include "ReplicatableVector.hpp"
00044 #include "ConstBoundaryCondition.hpp"
00045 #include "HeartEventHandler.hpp"
00046 #include "PetscMatTools.hpp"
00047 
00048 
00049 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00050 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::BoundaryConditionsContainer(bool deleteConditions)
00051             : AbstractBoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>(deleteConditions)
00052 {
00053     mLoadedFromArchive = false;
00054 
00055     for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00056     {
00057         mpNeumannMap[index_of_unknown] = new std::map< const BoundaryElement<ELEMENT_DIM-1, SPACE_DIM> *, const AbstractBoundaryCondition<SPACE_DIM>*>;
00058 
00059         mAnyNonZeroNeumannConditionsForUnknown[index_of_unknown] = false;
00060         mLastNeumannCondition[index_of_unknown] = mpNeumannMap[index_of_unknown]->begin();
00061 
00062         mpPeriodicBcMap[index_of_unknown] = new std::map< const Node<SPACE_DIM> *, const Node<SPACE_DIM> * >;
00063     }
00064 
00065     // This zero boundary condition is only used in AddNeumannBoundaryCondition
00066     mpZeroBoundaryCondition = new ConstBoundaryCondition<SPACE_DIM>(0.0);
00067 }
00068 
00069 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00070 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::~BoundaryConditionsContainer()
00071 {
00072     // Keep track of what boundary condition objects we've deleted
00073     std::set<const AbstractBoundaryCondition<SPACE_DIM>*> deleted_conditions;
00074     for (unsigned i=0; i<PROBLEM_DIM; i++)
00075     {
00076         NeumannMapIterator neumann_iterator = mpNeumannMap[i]->begin();
00077         while (neumann_iterator != mpNeumannMap[i]->end() )
00078         {
00079 
00080             if (deleted_conditions.count(neumann_iterator->second) == 0)
00081             {
00082                 deleted_conditions.insert(neumann_iterator->second);
00083                 //Leave the zero boundary condition until last
00084                 if (neumann_iterator->second != mpZeroBoundaryCondition)
00085                 {
00086                     if (this->mDeleteConditions)
00087                     {
00088                         delete neumann_iterator->second;
00089                     }
00090                 }
00091             }
00092             neumann_iterator++;
00093         }
00094         delete(mpNeumannMap[i]);
00095         delete(mpPeriodicBcMap[i]);
00096     }
00097 
00098     delete mpZeroBoundaryCondition;
00099 
00100     if (this->mDeleteConditions)
00101     {
00102         this->DeleteDirichletBoundaryConditions(deleted_conditions);
00103     }
00104 }
00105 
00106 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00107 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::AddDirichletBoundaryCondition(const Node<SPACE_DIM>* pBoundaryNode,
00108                                         const AbstractBoundaryCondition<SPACE_DIM>* pBoundaryCondition,
00109                                         unsigned indexOfUnknown,
00110                                         bool checkIfBoundaryNode)
00111 {
00112     assert(indexOfUnknown < PROBLEM_DIM);
00113     if (checkIfBoundaryNode)
00114     {
00115         assert(pBoundaryNode->IsBoundaryNode());
00116     }
00117 
00118     (*(this->mpDirichletMap[indexOfUnknown]))[pBoundaryNode] = pBoundaryCondition;
00119 }
00120 
00121 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00122 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::AddPeriodicBoundaryCondition(const Node<SPACE_DIM>* pNode1,
00123                                                                                                   const Node<SPACE_DIM>* pNode2)
00124 {
00125     assert(pNode1->IsBoundaryNode());
00126     assert(pNode2->IsBoundaryNode());
00127 
00128     // will assume the periodic BC is to be applied to ALL unknowns, can't really imagine a
00129     // situation where this isn't going to be true. If necessary can easily change this method
00130     // to take in the index of the unknown
00131     for(unsigned i=0; i<PROBLEM_DIM; i++)
00132     {
00133         (*(this->mpPeriodicBcMap[i]))[pNode1] = pNode2;
00134     }
00135 }
00136 
00137 
00138 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00139 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::AddNeumannBoundaryCondition( const BoundaryElement<ELEMENT_DIM-1, SPACE_DIM> * pBoundaryElement,
00140                                       const AbstractBoundaryCondition<SPACE_DIM> * pBoundaryCondition,
00141                                       unsigned indexOfUnknown)
00142 {
00143     assert(indexOfUnknown < PROBLEM_DIM);
00144 
00145     /*
00146      * If this condition is constant, we can test whether it is zero.
00147      * Otherwise we assume that this could be a non-zero boundary condition.
00148      */
00149     const ConstBoundaryCondition<SPACE_DIM>* p_const_cond = dynamic_cast<const ConstBoundaryCondition<SPACE_DIM>*>(pBoundaryCondition);
00150     if (p_const_cond)
00151     {
00152         if (p_const_cond->GetValue(pBoundaryElement->GetNode(0)->GetPoint()) != 0.0)
00153         {
00154             mAnyNonZeroNeumannConditionsForUnknown[indexOfUnknown] = true;
00155         }
00156     }
00157     else
00158     {
00159         mAnyNonZeroNeumannConditionsForUnknown[indexOfUnknown] = true;
00160     }
00161 
00162     for (unsigned unknown=0; unknown<PROBLEM_DIM; unknown++)
00163     {
00164         if (unknown==indexOfUnknown)
00165         {
00166             (*(mpNeumannMap[indexOfUnknown]))[pBoundaryElement] = pBoundaryCondition;
00167         }
00168         else
00169         {
00170             // If can't find pBoundaryElement in map[unknown]
00171             if ( mpNeumannMap[unknown]->find(pBoundaryElement)==mpNeumannMap[unknown]->end() )
00172             {
00173                 // Add zero bc to other unknowns (so all maps are in sync)
00174                 (*(mpNeumannMap[unknown]))[pBoundaryElement] = mpZeroBoundaryCondition;
00175             }
00176         }
00177     }
00178 }
00179 
00180 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00181 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::DefineZeroDirichletOnMeshBoundary(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00182                                            unsigned indexOfUnknown)
00183 {
00184     this->DefineConstantDirichletOnMeshBoundary(pMesh, 0.0, indexOfUnknown);
00185 }
00186 
00187 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00188 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::DefineConstantDirichletOnMeshBoundary(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00189                                                double value,
00190                                                unsigned indexOfUnknown)
00191 {
00192     assert(indexOfUnknown < PROBLEM_DIM);
00193 
00194     // In applying a condition to the boundary, we need to be sure that the boundary exists
00195     assert(pMesh->GetNumBoundaryNodes() > 0);
00196 
00197     ConstBoundaryCondition<SPACE_DIM>* p_boundary_condition = new ConstBoundaryCondition<SPACE_DIM>(value);
00198 
00199     typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::BoundaryNodeIterator iter;
00200     iter = pMesh->GetBoundaryNodeIteratorBegin();
00201     while (iter != pMesh->GetBoundaryNodeIteratorEnd())
00202     {
00203         AddDirichletBoundaryCondition(*iter, p_boundary_condition, indexOfUnknown);
00204         iter++;
00205     }
00206 }
00207 
00208 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00209 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::DefineZeroNeumannOnMeshBoundary(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00210                                          unsigned indexOfUnknown)
00211 {
00212     assert(indexOfUnknown < PROBLEM_DIM);
00213 
00214     // In applying a condition to the boundary, we need to be sure that the boundary exists
00215     assert(pMesh->GetNumBoundaryElements() > 0);
00216     ConstBoundaryCondition<SPACE_DIM>* p_zero_boundary_condition = new ConstBoundaryCondition<SPACE_DIM>( 0.0 );
00217 
00218     typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::BoundaryElementIterator iter;
00219     iter = pMesh->GetBoundaryElementIteratorBegin();
00220     while (iter != pMesh->GetBoundaryElementIteratorEnd())
00221     {
00222         AddNeumannBoundaryCondition(*iter, p_zero_boundary_condition, indexOfUnknown);
00223         iter++;
00224     }
00225 
00226     mAnyNonZeroNeumannConditionsForUnknown[indexOfUnknown] = false;
00227 }
00228 
00259 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00260 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyDirichletToLinearProblem(
00261         LinearSystem& rLinearSystem,
00262         bool applyToMatrix,
00263         bool applyToRhsVector)
00264 {
00265     HeartEventHandler::BeginEvent(HeartEventHandler::DIRICHLET_BCS);
00266 
00267     if (applyToMatrix)
00268     {
00269         if (!this->HasDirichletBoundaryConditions())
00270         {
00271             // Short-circuit the replication if there are no conditions
00272             HeartEventHandler::EndEvent(HeartEventHandler::DIRICHLET_BCS);
00273             return;
00274         }
00275 
00276         bool matrix_is_symmetric = rLinearSystem.IsMatrixSymmetric();
00277 
00278         if (matrix_is_symmetric)
00279         {
00280             /*
00281              * Modifications to the RHS are stored in the Dirichlet boundary
00282              * conditions vector. This is done so that they can be reapplied
00283              * at each time step.
00284              * Make a new vector to store the Dirichlet offsets in.
00285              */
00286             Vec& r_bcs_vec = rLinearSystem.rGetDirichletBoundaryConditionsVector();
00287             if (!r_bcs_vec)
00288             {
00289                 VecDuplicate(rLinearSystem.rGetRhsVector(), &r_bcs_vec);
00290             }
00291             PetscVecTools::Zero(r_bcs_vec);
00292             /*
00293              * If the matrix is symmetric, calls to GetMatrixRowDistributed()
00294              * require the matrix to be in assembled state. Otherwise we can
00295              * defer it.
00296              */
00297             rLinearSystem.AssembleFinalLinearSystem();
00298         }
00299 
00300         // Work out where we're setting Dirichlet boundary conditions *everywhere*, not just those locally known
00301         ReplicatableVector dirichlet_conditions(rLinearSystem.GetSize());
00302         unsigned lo, hi;
00303         {
00304             PetscInt lo_s, hi_s;
00305             rLinearSystem.GetOwnershipRange(lo_s, hi_s);
00306             lo = lo_s; hi = hi_s;
00307         }
00308         // Initialise all local entries to DBL_MAX, i.e. don't know if there's a condition
00309         for (unsigned i=lo; i<hi; i++)
00310         {
00311             dirichlet_conditions[i] = DBL_MAX;
00312         }
00313         // Now fill in the ones we know
00314         for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00315         {
00316             this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00317 
00318             while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00319             {
00320                 unsigned node_index = this->mDirichIterator->first->GetIndex();
00321                 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
00322                 assert(value != DBL_MAX);
00323 
00324                 unsigned row = PROBLEM_DIM*node_index + index_of_unknown;
00325                 dirichlet_conditions[row] = value;
00326 
00327                 this->mDirichIterator++;
00328             }
00329         }
00330 
00331         // And replicate
00332         dirichlet_conditions.Replicate(lo, hi);
00333 
00334         // Which rows have conditions?
00335         std::vector<unsigned> rows_to_zero;
00336         for (unsigned i=0; i<dirichlet_conditions.GetSize(); i++)
00337         {
00338             if (dirichlet_conditions[i] != DBL_MAX)
00339             {
00340                 rows_to_zero.push_back(i);
00341             }
00342         }
00343 
00344         if (matrix_is_symmetric)
00345         {
00346             // Modify the matrix columns
00347             for (unsigned i=0; i<rows_to_zero.size(); i++)
00348             {
00349                 unsigned col = rows_to_zero[i];
00350                 double minus_value = -dirichlet_conditions[col];
00351 
00352                 /*
00353                  * Get a vector which will store the column of the matrix (column d,
00354                  * where d is the index of the row (and column) to be altered for the
00355                  * boundary condition. Since the matrix is symmetric when get row
00356                  * number "col" and treat it as a column. PETSc uses compressed row
00357                  * format and therefore getting rows is far more efficient than getting
00358                  * columns.
00359                  */
00360                 Vec matrix_col = rLinearSystem.GetMatrixRowDistributed(col);
00361 
00362                 // Zero the correct entry of the column
00363                 PetscVecTools::SetElement(matrix_col, col, 0.0);
00364 
00365                 /*
00366                  * Set up the RHS Dirichlet boundary conditions vector.
00367                  * Assuming one boundary at the zeroth node (x_0 = value), this is equal to
00368                  *   -value*[0 a_21 a_31 .. a_N1]
00369                  * and will be added to the RHS.
00370                  */
00371                 PetscVecTools::AddScaledVector(rLinearSystem.rGetDirichletBoundaryConditionsVector(), matrix_col, minus_value);
00372                 PetscTools::Destroy(matrix_col);
00373             }
00374         }
00375 
00376         /*
00377          * Now zero the appropriate rows and columns of the matrix. If the matrix
00378          * is symmetric we apply the boundary conditions in a way the symmetry isn't
00379          * lost (rows and columns). If not only the row is zeroed.
00380          */
00381         if (matrix_is_symmetric)
00382         {
00383             rLinearSystem.ZeroMatrixRowsAndColumnsWithValueOnDiagonal(rows_to_zero, 1.0);
00384         }
00385         else
00386         {
00387             rLinearSystem.ZeroMatrixRowsWithValueOnDiagonal(rows_to_zero, 1.0);
00388         }
00389     }
00390 
00391     if (applyToRhsVector)
00392     {
00393         // Apply the RHS boundary conditions modification if required.
00394         if (rLinearSystem.rGetDirichletBoundaryConditionsVector())
00395         {
00396             PetscVecTools::AddScaledVector(rLinearSystem.rGetRhsVector(), rLinearSystem.rGetDirichletBoundaryConditionsVector(), 1.0);
00397         }
00398 
00399         /*
00400          * Apply the actual boundary condition to the RHS, note this must be
00401          * done after the modification to the RHS vector.
00402          */
00403         for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00404         {
00405             this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00406 
00407             while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00408             {
00409                 unsigned node_index = this->mDirichIterator->first->GetIndex();
00410                 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
00411 
00412                 unsigned row = PROBLEM_DIM*node_index + index_of_unknown;
00413 
00414                 rLinearSystem.SetRhsVectorElement(row, value);
00415 
00416                 this->mDirichIterator++;
00417             }
00418         }
00419     }
00420 
00421     HeartEventHandler::EndEvent(HeartEventHandler::DIRICHLET_BCS);
00422 }
00423 
00424 
00425 
00426 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00427 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyPeriodicBcsToLinearProblem(LinearSystem& rLinearSystem,
00428                                                                                                      bool applyToMatrix,
00429                                                                                                      bool applyToRhsVector)
00430 {
00431     bool has_periodic_bcs = false;
00432     for (unsigned i=0; i<PROBLEM_DIM; i++)
00433     {
00434         if (!mpPeriodicBcMap[i]->empty())
00435         {
00436             has_periodic_bcs = true;
00437             break;
00438         }
00439     }
00440 
00441     if(!has_periodic_bcs)
00442     {
00443         return;
00444     }
00445 
00446     if(applyToMatrix)
00447     {
00448         std::vector<unsigned> rows_to_zero;
00449         for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00450         {
00451             for(typename std::map< const Node<SPACE_DIM> *, const Node<SPACE_DIM> * >::const_iterator iter = mpPeriodicBcMap[index_of_unknown]->begin();
00452                 iter != mpPeriodicBcMap[index_of_unknown]->end();
00453                 ++iter)
00454             {
00455                 unsigned node_index_1 = iter->first->GetIndex();
00456                 unsigned row_index_1 = PROBLEM_DIM*node_index_1 + index_of_unknown;
00457                 rows_to_zero.push_back(row_index_1);
00458             }
00459         }
00460 
00461         rLinearSystem.ZeroMatrixRowsWithValueOnDiagonal(rows_to_zero, 1.0);
00462 
00463         for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00464         {
00465             for(typename std::map< const Node<SPACE_DIM> *, const Node<SPACE_DIM> * >::const_iterator iter = mpPeriodicBcMap[index_of_unknown]->begin();
00466                 iter != mpPeriodicBcMap[index_of_unknown]->end();
00467                 ++iter)
00468             {
00469                 unsigned node_index_1 = iter->first->GetIndex();
00470                 unsigned node_index_2 = iter->second->GetIndex();
00471 
00472                 unsigned mat_index1 = PROBLEM_DIM*node_index_1 + index_of_unknown;
00473                 unsigned mat_index2 = PROBLEM_DIM*node_index_2 + index_of_unknown;
00474                 PetscMatTools::SetElement(rLinearSystem.rGetLhsMatrix(), mat_index1, mat_index2, -1.0);
00475             }
00476         }
00477     }
00478 
00479     if(applyToRhsVector)
00480     {
00481         for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00482         {
00483             for(typename std::map< const Node<SPACE_DIM> *, const Node<SPACE_DIM> * >::const_iterator iter = mpPeriodicBcMap[index_of_unknown]->begin();
00484                 iter != mpPeriodicBcMap[index_of_unknown]->end();
00485                 ++iter)
00486             {
00487                 unsigned node_index = iter->first->GetIndex();
00488                 unsigned row_index = PROBLEM_DIM*node_index + index_of_unknown;
00489                 rLinearSystem.SetRhsVectorElement(row_index, 0.0);
00490             }
00491         }
00492     }
00493 }
00494 
00495 
00496 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00497 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyDirichletToNonlinearResidual(
00498         const Vec currentSolution,
00499         Vec residual,
00500         DistributedVectorFactory& rFactory)
00501 {
00502     for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00503     {
00504         this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00505 
00506         DistributedVector solution_distributed = rFactory.CreateDistributedVector(currentSolution);
00507         DistributedVector residual_distributed = rFactory.CreateDistributedVector(residual);
00508 
00509 
00510         while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00511         {
00512             DistributedVector::Stripe solution_stripe(solution_distributed, index_of_unknown);
00513             DistributedVector::Stripe residual_stripe(residual_distributed, index_of_unknown);
00514 
00515             unsigned node_index = this->mDirichIterator->first->GetIndex();
00516 
00517             double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
00518 
00519             if (solution_distributed.IsGlobalIndexLocal(node_index))
00520             {
00521                 residual_stripe[node_index]=solution_stripe[node_index] - value;
00522             }
00523             this->mDirichIterator++;
00524         }
00525         solution_distributed.Restore();
00526         residual_distributed.Restore();
00527     }
00528 }
00529 
00530 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00531 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyDirichletToNonlinearJacobian(Mat jacobian)
00532 {
00533     unsigned num_boundary_conditions = 0;
00534     for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00535     {
00536         num_boundary_conditions += this->mpDirichletMap[index_of_unknown]->size();
00537     }
00538 
00539     std::vector<unsigned> rows_to_zero(num_boundary_conditions);
00540 
00541     unsigned counter=0;
00542     for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00543     {
00544         this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00545 
00546         while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00547         {
00548             unsigned node_index = this->mDirichIterator->first->GetIndex();
00549             rows_to_zero[counter++] = PROBLEM_DIM*node_index + index_of_unknown;
00550             this->mDirichIterator++;
00551         }
00552     }
00553     PetscMatTools::Finalise(jacobian);
00554     PetscMatTools::ZeroRowsWithValueOnDiagonal(jacobian, rows_to_zero, 1.0);
00555 }
00556 
00557 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00558 bool BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::Validate(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh)
00559 {
00560     bool valid = true;
00561 
00562     for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00563     {
00564         // Iterate over surface elements
00565         typename AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>::BoundaryElementIterator elt_iter
00566         = pMesh->GetBoundaryElementIteratorBegin();
00567         while (valid && elt_iter != pMesh->GetBoundaryElementIteratorEnd())
00568         {
00569             if (!HasNeumannBoundaryCondition(*elt_iter, index_of_unknown))
00570             {
00571                 // Check for Dirichlet conditions on this element's nodes
00572                 for (unsigned i=0; i<(*elt_iter)->GetNumNodes(); i++)
00573                 {
00574                     if (!this->HasDirichletBoundaryCondition((*elt_iter)->GetNode(i)))
00575                     {
00576                         valid = false;
00577                     }
00578                 }
00579             }
00580             elt_iter++;
00581         }
00582     }
00583     return valid;
00584 }
00585 
00586 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00587 double BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::GetNeumannBCValue(const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>* pSurfaceElement,
00588                              const ChastePoint<SPACE_DIM>& rX,
00589                              unsigned indexOfUnknown)
00590 {
00591     assert(indexOfUnknown < PROBLEM_DIM);
00592 
00593     // Did we see this condition on the last search we did?
00594     if (mLastNeumannCondition[indexOfUnknown] == mpNeumannMap[indexOfUnknown]->end() ||
00595         mLastNeumannCondition[indexOfUnknown]->first != pSurfaceElement)
00596     {
00597         mLastNeumannCondition[indexOfUnknown] = mpNeumannMap[indexOfUnknown]->find(pSurfaceElement);
00598     }
00599     if (mLastNeumannCondition[indexOfUnknown] == mpNeumannMap[indexOfUnknown]->end())
00600     {
00601         // No Neumann condition is equivalent to a zero Neumann condition
00602         return 0.0;
00603     }
00604     else
00605     {
00606         return mLastNeumannCondition[indexOfUnknown]->second->GetValue(rX);
00607     }
00608 }
00609 
00610 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00611 bool BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::HasNeumannBoundaryCondition(const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>* pSurfaceElement,
00612                                      unsigned indexOfUnknown)
00613 {
00614     assert(indexOfUnknown < PROBLEM_DIM);
00615 
00616     mLastNeumannCondition[indexOfUnknown] = mpNeumannMap[indexOfUnknown]->find(pSurfaceElement);
00617 
00618     return (mLastNeumannCondition[indexOfUnknown] != mpNeumannMap[indexOfUnknown]->end());
00619 }
00620 
00621 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00622 bool BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::AnyNonZeroNeumannConditions()
00623 {
00624     bool ret = false;
00625     for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00626     {
00627         if (mAnyNonZeroNeumannConditionsForUnknown[index_of_unknown] == true)
00628         {
00629             ret = true;
00630         }
00631     }
00632     return ret;
00633 }
00634 
00635 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00636 typename BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::NeumannMapIterator BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::BeginNeumann()
00637 {
00638     // [0] is ok as all maps will be in sync due to the way ApplyNeumannBoundaryCondition works
00639     return mpNeumannMap[0]->begin();
00640 }
00641 
00642 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00643 typename BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::NeumannMapIterator BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::EndNeumann()
00644 {
00645     // [0] is ok as all maps will be in sync due to the way ApplyNeumannBoundaryCondition works
00646     return mpNeumannMap[0]->end();
00647 }
00648 
00649 #endif // _BOUNDARYCONDITIONSCONTAINERIMPLEMENTATION_HPP_