00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
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
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
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
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
00129
00130
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
00147
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
00171 if ( mpNeumannMap[unknown]->find(pBoundaryElement)==mpNeumannMap[unknown]->end() )
00172 {
00173
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 assert(PetscTools::ReplicateBool( pMesh->GetNumBoundaryNodes() > 0 ) );
00195
00196 ConstBoundaryCondition<SPACE_DIM>* p_boundary_condition = new ConstBoundaryCondition<SPACE_DIM>(value);
00197
00198 typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::BoundaryNodeIterator iter;
00199 iter = pMesh->GetBoundaryNodeIteratorBegin();
00200 while (iter != pMesh->GetBoundaryNodeIteratorEnd())
00201 {
00202 AddDirichletBoundaryCondition(*iter, p_boundary_condition, indexOfUnknown);
00203 iter++;
00204 }
00205 }
00206
00207 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00208 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::DefineZeroNeumannOnMeshBoundary(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00209 unsigned indexOfUnknown)
00210 {
00211 assert(indexOfUnknown < PROBLEM_DIM);
00212
00213
00214 assert(pMesh->GetNumBoundaryElements() > 0);
00215 ConstBoundaryCondition<SPACE_DIM>* p_zero_boundary_condition = new ConstBoundaryCondition<SPACE_DIM>( 0.0 );
00216
00217 typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::BoundaryElementIterator iter;
00218 iter = pMesh->GetBoundaryElementIteratorBegin();
00219 while (iter != pMesh->GetBoundaryElementIteratorEnd())
00220 {
00221 AddNeumannBoundaryCondition(*iter, p_zero_boundary_condition, indexOfUnknown);
00222 iter++;
00223 }
00224
00225 mAnyNonZeroNeumannConditionsForUnknown[indexOfUnknown] = false;
00226 }
00227
00258 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00259 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyDirichletToLinearProblem(
00260 LinearSystem& rLinearSystem,
00261 bool applyToMatrix,
00262 bool applyToRhsVector)
00263 {
00264 HeartEventHandler::BeginEvent(HeartEventHandler::DIRICHLET_BCS);
00265
00266 if (applyToMatrix)
00267 {
00268 if (!this->HasDirichletBoundaryConditions())
00269 {
00270
00271 HeartEventHandler::EndEvent(HeartEventHandler::DIRICHLET_BCS);
00272 return;
00273 }
00274
00275 bool matrix_is_symmetric = rLinearSystem.IsMatrixSymmetric();
00276
00277 if (matrix_is_symmetric)
00278 {
00279
00280
00281
00282
00283
00284
00285 Vec& r_bcs_vec = rLinearSystem.rGetDirichletBoundaryConditionsVector();
00286 if (!r_bcs_vec)
00287 {
00288 VecDuplicate(rLinearSystem.rGetRhsVector(), &r_bcs_vec);
00289 }
00290 PetscVecTools::Zero(r_bcs_vec);
00291
00292
00293
00294
00295
00296 rLinearSystem.AssembleFinalLinearSystem();
00297 }
00298
00299
00300 ReplicatableVector dirichlet_conditions(rLinearSystem.GetSize());
00301 unsigned lo, hi;
00302 {
00303 PetscInt lo_s, hi_s;
00304 rLinearSystem.GetOwnershipRange(lo_s, hi_s);
00305 lo = lo_s; hi = hi_s;
00306 }
00307
00308 for (unsigned i=lo; i<hi; i++)
00309 {
00310 dirichlet_conditions[i] = DBL_MAX;
00311 }
00312
00313 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00314 {
00315 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00316
00317 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00318 {
00319 unsigned node_index = this->mDirichIterator->first->GetIndex();
00320 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
00321 assert(value != DBL_MAX);
00322
00323 unsigned row = PROBLEM_DIM*node_index + index_of_unknown;
00324 dirichlet_conditions[row] = value;
00325
00326 this->mDirichIterator++;
00327 }
00328 }
00329
00330
00331 dirichlet_conditions.Replicate(lo, hi);
00332
00333
00334 std::vector<unsigned> rows_to_zero;
00335 for (unsigned i=0; i<dirichlet_conditions.GetSize(); i++)
00336 {
00337 if (dirichlet_conditions[i] != DBL_MAX)
00338 {
00339 rows_to_zero.push_back(i);
00340 }
00341 }
00342
00343 if (matrix_is_symmetric)
00344 {
00345
00346 for (unsigned i=0; i<rows_to_zero.size(); i++)
00347 {
00348 unsigned col = rows_to_zero[i];
00349 double minus_value = -dirichlet_conditions[col];
00350
00351
00352
00353
00354
00355
00356
00357
00358
00359 Vec matrix_col = rLinearSystem.GetMatrixRowDistributed(col);
00360
00361
00362 PetscVecTools::SetElement(matrix_col, col, 0.0);
00363
00364
00365
00366
00367
00368
00369
00370 PetscVecTools::AddScaledVector(rLinearSystem.rGetDirichletBoundaryConditionsVector(), matrix_col, minus_value);
00371 PetscTools::Destroy(matrix_col);
00372 }
00373 }
00374
00375
00376
00377
00378
00379
00380 if (matrix_is_symmetric)
00381 {
00382 rLinearSystem.ZeroMatrixRowsAndColumnsWithValueOnDiagonal(rows_to_zero, 1.0);
00383 }
00384 else
00385 {
00386 rLinearSystem.ZeroMatrixRowsWithValueOnDiagonal(rows_to_zero, 1.0);
00387 }
00388 }
00389
00390 if (applyToRhsVector)
00391 {
00392
00393 if (rLinearSystem.rGetDirichletBoundaryConditionsVector())
00394 {
00395 PetscVecTools::AddScaledVector(rLinearSystem.rGetRhsVector(), rLinearSystem.rGetDirichletBoundaryConditionsVector(), 1.0);
00396 }
00397
00398
00399
00400
00401
00402 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00403 {
00404 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00405
00406 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00407 {
00408 unsigned node_index = this->mDirichIterator->first->GetIndex();
00409 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
00410
00411 unsigned row = PROBLEM_DIM*node_index + index_of_unknown;
00412
00413 rLinearSystem.SetRhsVectorElement(row, value);
00414
00415 this->mDirichIterator++;
00416 }
00417 }
00418 }
00419
00420 HeartEventHandler::EndEvent(HeartEventHandler::DIRICHLET_BCS);
00421 }
00422
00423
00424
00425 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00426 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyPeriodicBcsToLinearProblem(LinearSystem& rLinearSystem,
00427 bool applyToMatrix,
00428 bool applyToRhsVector)
00429 {
00430 bool has_periodic_bcs = false;
00431 for (unsigned i=0; i<PROBLEM_DIM; i++)
00432 {
00433 if (!mpPeriodicBcMap[i]->empty())
00434 {
00435 has_periodic_bcs = true;
00436 break;
00437 }
00438 }
00439
00440 if(!has_periodic_bcs)
00441 {
00442 return;
00443 }
00444
00445 if(applyToMatrix)
00446 {
00447 std::vector<unsigned> rows_to_zero;
00448 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00449 {
00450 for(typename std::map< const Node<SPACE_DIM> *, const Node<SPACE_DIM> * >::const_iterator iter = mpPeriodicBcMap[index_of_unknown]->begin();
00451 iter != mpPeriodicBcMap[index_of_unknown]->end();
00452 ++iter)
00453 {
00454 unsigned node_index_1 = iter->first->GetIndex();
00455 unsigned row_index_1 = PROBLEM_DIM*node_index_1 + index_of_unknown;
00456 rows_to_zero.push_back(row_index_1);
00457 }
00458 }
00459
00460 rLinearSystem.ZeroMatrixRowsWithValueOnDiagonal(rows_to_zero, 1.0);
00461
00462 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00463 {
00464 for(typename std::map< const Node<SPACE_DIM> *, const Node<SPACE_DIM> * >::const_iterator iter = mpPeriodicBcMap[index_of_unknown]->begin();
00465 iter != mpPeriodicBcMap[index_of_unknown]->end();
00466 ++iter)
00467 {
00468 unsigned node_index_1 = iter->first->GetIndex();
00469 unsigned node_index_2 = iter->second->GetIndex();
00470
00471 unsigned mat_index1 = PROBLEM_DIM*node_index_1 + index_of_unknown;
00472 unsigned mat_index2 = PROBLEM_DIM*node_index_2 + index_of_unknown;
00473 PetscMatTools::SetElement(rLinearSystem.rGetLhsMatrix(), mat_index1, mat_index2, -1.0);
00474 }
00475 }
00476 }
00477
00478 if(applyToRhsVector)
00479 {
00480 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00481 {
00482 for(typename std::map< const Node<SPACE_DIM> *, const Node<SPACE_DIM> * >::const_iterator iter = mpPeriodicBcMap[index_of_unknown]->begin();
00483 iter != mpPeriodicBcMap[index_of_unknown]->end();
00484 ++iter)
00485 {
00486 unsigned node_index = iter->first->GetIndex();
00487 unsigned row_index = PROBLEM_DIM*node_index + index_of_unknown;
00488 rLinearSystem.SetRhsVectorElement(row_index, 0.0);
00489 }
00490 }
00491 }
00492 }
00493
00494
00495 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00496 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyDirichletToNonlinearResidual(
00497 const Vec currentSolution,
00498 Vec residual,
00499 DistributedVectorFactory& rFactory)
00500 {
00501 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00502 {
00503 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00504
00505 DistributedVector solution_distributed = rFactory.CreateDistributedVector(currentSolution);
00506 DistributedVector residual_distributed = rFactory.CreateDistributedVector(residual);
00507
00508
00509 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00510 {
00511 DistributedVector::Stripe solution_stripe(solution_distributed, index_of_unknown);
00512 DistributedVector::Stripe residual_stripe(residual_distributed, index_of_unknown);
00513
00514 unsigned node_index = this->mDirichIterator->first->GetIndex();
00515
00516 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
00517
00518 if (solution_distributed.IsGlobalIndexLocal(node_index))
00519 {
00520 residual_stripe[node_index]=solution_stripe[node_index] - value;
00521 }
00522 this->mDirichIterator++;
00523 }
00524 solution_distributed.Restore();
00525 residual_distributed.Restore();
00526 }
00527 }
00528
00529 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00530 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyDirichletToNonlinearJacobian(Mat jacobian)
00531 {
00532 unsigned num_boundary_conditions = 0;
00533 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00534 {
00535 num_boundary_conditions += this->mpDirichletMap[index_of_unknown]->size();
00536 }
00537
00538 std::vector<unsigned> rows_to_zero(num_boundary_conditions);
00539
00540 unsigned counter=0;
00541 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00542 {
00543 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00544
00545 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00546 {
00547 unsigned node_index = this->mDirichIterator->first->GetIndex();
00548 rows_to_zero[counter++] = PROBLEM_DIM*node_index + index_of_unknown;
00549 this->mDirichIterator++;
00550 }
00551 }
00552 PetscMatTools::Finalise(jacobian);
00553 PetscMatTools::ZeroRowsWithValueOnDiagonal(jacobian, rows_to_zero, 1.0);
00554 }
00555
00556 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00557 bool BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::Validate(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh)
00558 {
00559 bool valid = true;
00560
00561 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00562 {
00563
00564 typename AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>::BoundaryElementIterator elt_iter
00565 = pMesh->GetBoundaryElementIteratorBegin();
00566 while (valid && elt_iter != pMesh->GetBoundaryElementIteratorEnd())
00567 {
00568 if (!HasNeumannBoundaryCondition(*elt_iter, index_of_unknown))
00569 {
00570
00571 for (unsigned i=0; i<(*elt_iter)->GetNumNodes(); i++)
00572 {
00573 if (!this->HasDirichletBoundaryCondition((*elt_iter)->GetNode(i)))
00574 {
00575 valid = false;
00576 }
00577 }
00578 }
00579 elt_iter++;
00580 }
00581 }
00582 return valid;
00583 }
00584
00585 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00586 double BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::GetNeumannBCValue(const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>* pSurfaceElement,
00587 const ChastePoint<SPACE_DIM>& rX,
00588 unsigned indexOfUnknown)
00589 {
00590 assert(indexOfUnknown < PROBLEM_DIM);
00591
00592
00593 if (mLastNeumannCondition[indexOfUnknown] == mpNeumannMap[indexOfUnknown]->end() ||
00594 mLastNeumannCondition[indexOfUnknown]->first != pSurfaceElement)
00595 {
00596 mLastNeumannCondition[indexOfUnknown] = mpNeumannMap[indexOfUnknown]->find(pSurfaceElement);
00597 }
00598 if (mLastNeumannCondition[indexOfUnknown] == mpNeumannMap[indexOfUnknown]->end())
00599 {
00600
00601 return 0.0;
00602 }
00603 else
00604 {
00605 return mLastNeumannCondition[indexOfUnknown]->second->GetValue(rX);
00606 }
00607 }
00608
00609 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00610 bool BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::HasNeumannBoundaryCondition(const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>* pSurfaceElement,
00611 unsigned indexOfUnknown)
00612 {
00613 assert(indexOfUnknown < PROBLEM_DIM);
00614
00615 mLastNeumannCondition[indexOfUnknown] = mpNeumannMap[indexOfUnknown]->find(pSurfaceElement);
00616
00617 return (mLastNeumannCondition[indexOfUnknown] != mpNeumannMap[indexOfUnknown]->end());
00618 }
00619
00620 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00621 bool BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::AnyNonZeroNeumannConditions()
00622 {
00623 bool ret = false;
00624 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00625 {
00626 if (mAnyNonZeroNeumannConditionsForUnknown[index_of_unknown] == true)
00627 {
00628 ret = true;
00629 }
00630 }
00631 return ret;
00632 }
00633
00634 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00635 typename BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::NeumannMapIterator BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::BeginNeumann()
00636 {
00637
00638 return mpNeumannMap[0]->begin();
00639 }
00640
00641 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00642 typename BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::NeumannMapIterator BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::EndNeumann()
00643 {
00644
00645 return mpNeumannMap[0]->end();
00646 }
00647
00648 #endif // _BOUNDARYCONDITIONSCONTAINERIMPLEMENTATION_HPP_