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 #ifndef _BOUNDARYCONDITIONSCONTAINERIMPLEMENTATION_HPP_
00030 #define _BOUNDARYCONDITIONSCONTAINERIMPLEMENTATION_HPP_
00031
00032 #include "PetscVecTools.hpp"
00033 #include "PetscMatTools.hpp"
00034 #include "BoundaryConditionsContainer.hpp"
00035 #include "DistributedVector.hpp"
00036 #include "ReplicatableVector.hpp"
00037 #include "ConstBoundaryCondition.hpp"
00038
00039 #include "HeartEventHandler.hpp"
00040
00041 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00042 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::BoundaryConditionsContainer(bool deleteConditions)
00043 : AbstractBoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>(deleteConditions)
00044 {
00045 mLoadedFromArchive = false;
00046
00047 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00048 {
00049 mpNeumannMap[index_of_unknown] = new std::map< const BoundaryElement<ELEMENT_DIM-1, SPACE_DIM> *, const AbstractBoundaryCondition<SPACE_DIM>*>;
00050
00051 mAnyNonZeroNeumannConditionsForUnknown[index_of_unknown] = false;
00052 mLastNeumannCondition[index_of_unknown] = mpNeumannMap[index_of_unknown]->begin();
00053 }
00054
00055
00056 mpZeroBoundaryCondition = new ConstBoundaryCondition<SPACE_DIM>(0.0);
00057 }
00058
00059 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00060 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::~BoundaryConditionsContainer()
00061 {
00062
00063 std::set<const AbstractBoundaryCondition<SPACE_DIM>*> deleted_conditions;
00064 for (unsigned i=0; i<PROBLEM_DIM; i++)
00065 {
00066 NeumannMapIterator neumann_iterator = mpNeumannMap[i]->begin();
00067 while (neumann_iterator != mpNeumannMap[i]->end() )
00068 {
00069
00070 if (deleted_conditions.count(neumann_iterator->second) == 0)
00071 {
00072 deleted_conditions.insert(neumann_iterator->second);
00073
00074 if (neumann_iterator->second != mpZeroBoundaryCondition)
00075 {
00076 if (this->mDeleteConditions)
00077 {
00078 delete neumann_iterator->second;
00079 }
00080 }
00081 }
00082 neumann_iterator++;
00083 }
00084 delete(mpNeumannMap[i]);
00085 }
00086
00087 delete mpZeroBoundaryCondition;
00088
00089 if (this->mDeleteConditions)
00090 {
00091 this->DeleteDirichletBoundaryConditions(deleted_conditions);
00092 }
00093 }
00094
00095 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00096 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::AddDirichletBoundaryCondition(const Node<SPACE_DIM>* pBoundaryNode,
00097 const AbstractBoundaryCondition<SPACE_DIM>* pBoundaryCondition,
00098 unsigned indexOfUnknown,
00099 bool checkIfBoundaryNode)
00100 {
00101 assert(indexOfUnknown < PROBLEM_DIM);
00102 if (checkIfBoundaryNode)
00103 {
00104 assert(pBoundaryNode->IsBoundaryNode());
00105 }
00106
00107 (*(this->mpDirichletMap[indexOfUnknown]))[pBoundaryNode] = pBoundaryCondition;
00108 }
00109
00110 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00111 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::AddNeumannBoundaryCondition( const BoundaryElement<ELEMENT_DIM-1, SPACE_DIM> * pBoundaryElement,
00112 const AbstractBoundaryCondition<SPACE_DIM> * pBoundaryCondition,
00113 unsigned indexOfUnknown)
00114 {
00115 assert(indexOfUnknown < PROBLEM_DIM);
00116
00117
00118
00119
00120
00121 const ConstBoundaryCondition<SPACE_DIM>* p_const_cond = dynamic_cast<const ConstBoundaryCondition<SPACE_DIM>*>(pBoundaryCondition);
00122 if (p_const_cond)
00123 {
00124 if (p_const_cond->GetValue(pBoundaryElement->GetNode(0)->GetPoint()) != 0.0)
00125 {
00126 mAnyNonZeroNeumannConditionsForUnknown[indexOfUnknown] = true;
00127 }
00128 }
00129 else
00130 {
00131 mAnyNonZeroNeumannConditionsForUnknown[indexOfUnknown] = true;
00132 }
00133
00134 for (unsigned unknown=0; unknown<PROBLEM_DIM; unknown++)
00135 {
00136 if (unknown==indexOfUnknown)
00137 {
00138 (*(mpNeumannMap[indexOfUnknown]))[pBoundaryElement] = pBoundaryCondition;
00139 }
00140 else
00141 {
00142
00143 if ( mpNeumannMap[unknown]->find(pBoundaryElement)==mpNeumannMap[unknown]->end() )
00144 {
00145
00146 (*(mpNeumannMap[unknown]))[pBoundaryElement] = mpZeroBoundaryCondition;
00147 }
00148 }
00149 }
00150 }
00151
00152 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00153 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::DefineZeroDirichletOnMeshBoundary(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00154 unsigned indexOfUnknown)
00155 {
00156 this->DefineConstantDirichletOnMeshBoundary(pMesh, 0.0, indexOfUnknown);
00157 }
00158
00159 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00160 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::DefineConstantDirichletOnMeshBoundary(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00161 double value,
00162 unsigned indexOfUnknown)
00163 {
00164 assert(indexOfUnknown < PROBLEM_DIM);
00165
00166
00167 assert(pMesh->GetNumBoundaryNodes() > 0);
00168
00169 ConstBoundaryCondition<SPACE_DIM>* p_boundary_condition = new ConstBoundaryCondition<SPACE_DIM>(value);
00170
00171 typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::BoundaryNodeIterator iter;
00172 iter = pMesh->GetBoundaryNodeIteratorBegin();
00173 while (iter != pMesh->GetBoundaryNodeIteratorEnd())
00174 {
00175 AddDirichletBoundaryCondition(*iter, p_boundary_condition, indexOfUnknown);
00176 iter++;
00177 }
00178 }
00179
00180 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00181 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::DefineZeroNeumannOnMeshBoundary(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00182 unsigned indexOfUnknown)
00183 {
00184 assert(indexOfUnknown < PROBLEM_DIM);
00185
00186
00187 assert(pMesh->GetNumBoundaryElements() > 0);
00188 ConstBoundaryCondition<SPACE_DIM>* p_zero_boundary_condition = new ConstBoundaryCondition<SPACE_DIM>( 0.0 );
00189
00190 typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::BoundaryElementIterator iter;
00191 iter = pMesh->GetBoundaryElementIteratorBegin();
00192 while (iter != pMesh->GetBoundaryElementIteratorEnd())
00193 {
00194 AddNeumannBoundaryCondition(*iter, p_zero_boundary_condition, indexOfUnknown);
00195 iter++;
00196 }
00197
00198 mAnyNonZeroNeumannConditionsForUnknown[indexOfUnknown] = false;
00199 }
00200
00231 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00232 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyDirichletToLinearProblem(
00233 LinearSystem& rLinearSystem,
00234 bool applyToMatrix,
00235 bool applyToRhsVector)
00236 {
00237 HeartEventHandler::BeginEvent(HeartEventHandler::DIRICHLET_BCS);
00238
00239 if (applyToMatrix)
00240 {
00241 if (!this->HasDirichletBoundaryConditions())
00242 {
00243
00244 HeartEventHandler::EndEvent(HeartEventHandler::DIRICHLET_BCS);
00245 return;
00246 }
00247
00248 bool matrix_is_symmetric = rLinearSystem.IsMatrixSymmetric();
00249
00250 if (matrix_is_symmetric)
00251 {
00252
00253
00254
00255
00256
00257
00258 Vec& r_bcs_vec = rLinearSystem.rGetDirichletBoundaryConditionsVector();
00259 if (!r_bcs_vec)
00260 {
00261 VecDuplicate(rLinearSystem.rGetRhsVector(), &r_bcs_vec);
00262 }
00263 PetscVecTools::Zero(r_bcs_vec);
00264
00265
00266
00267
00268
00269 rLinearSystem.AssembleFinalLinearSystem();
00270 }
00271
00272
00273 ReplicatableVector dirichlet_conditions(rLinearSystem.GetSize());
00274 unsigned lo, hi;
00275 {
00276 PetscInt lo_s, hi_s;
00277 rLinearSystem.GetOwnershipRange(lo_s, hi_s);
00278 lo = lo_s; hi = hi_s;
00279 }
00280
00281 for (unsigned i=lo; i<hi; i++)
00282 {
00283 dirichlet_conditions[i] = DBL_MAX;
00284 }
00285
00286 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00287 {
00288 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00289
00290 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00291 {
00292 unsigned node_index = this->mDirichIterator->first->GetIndex();
00293 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
00294 assert(value != DBL_MAX);
00295
00296 unsigned row = PROBLEM_DIM*node_index + index_of_unknown;
00297 dirichlet_conditions[row] = value;
00298
00299 this->mDirichIterator++;
00300 }
00301 }
00302
00303
00304 dirichlet_conditions.Replicate(lo, hi);
00305
00306
00307 std::vector<unsigned> rows_to_zero;
00308 for (unsigned i=0; i<dirichlet_conditions.GetSize(); i++)
00309 {
00310 if (dirichlet_conditions[i] != DBL_MAX)
00311 {
00312 rows_to_zero.push_back(i);
00313 }
00314 }
00315
00316 if (matrix_is_symmetric)
00317 {
00318
00319 for (unsigned i=0; i<rows_to_zero.size(); i++)
00320 {
00321 unsigned col = rows_to_zero[i];
00322 double minus_value = -dirichlet_conditions[col];
00323
00324
00325
00326
00327
00328
00329
00330
00331
00332 Vec matrix_col = rLinearSystem.GetMatrixRowDistributed(col);
00333
00334
00335 PetscVecTools::SetElement(matrix_col, col, 0.0);
00336
00337
00338
00339
00340
00341
00342
00343 PetscVecTools::AddScaledVector(rLinearSystem.rGetDirichletBoundaryConditionsVector(), matrix_col, minus_value);
00344 VecDestroy(matrix_col);
00345 }
00346 }
00347
00348
00349
00350
00351
00352
00353 if (matrix_is_symmetric)
00354 {
00355 rLinearSystem.ZeroMatrixRowsAndColumnsWithValueOnDiagonal(rows_to_zero, 1.0);
00356 }
00357 else
00358 {
00359 rLinearSystem.ZeroMatrixRowsWithValueOnDiagonal(rows_to_zero, 1.0);
00360 }
00361 }
00362
00363 if (applyToRhsVector)
00364 {
00365
00366 if (rLinearSystem.rGetDirichletBoundaryConditionsVector())
00367 {
00368 PetscVecTools::AddScaledVector(rLinearSystem.rGetRhsVector(), rLinearSystem.rGetDirichletBoundaryConditionsVector(), 1.0);
00369 }
00370
00371
00372
00373
00374
00375 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00376 {
00377 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00378
00379 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00380 {
00381 unsigned node_index = this->mDirichIterator->first->GetIndex();
00382 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
00383
00384 unsigned row = PROBLEM_DIM*node_index + index_of_unknown;
00385
00386 rLinearSystem.SetRhsVectorElement(row, value);
00387
00388 this->mDirichIterator++;
00389 }
00390 }
00391 }
00392
00393 HeartEventHandler::EndEvent(HeartEventHandler::DIRICHLET_BCS);
00394 }
00395
00396 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00397 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyDirichletToNonlinearResidual(
00398 const Vec currentSolution,
00399 Vec residual,
00400 DistributedVectorFactory& rFactory)
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 DistributedVector solution_distributed = rFactory.CreateDistributedVector(currentSolution);
00407 DistributedVector residual_distributed = rFactory.CreateDistributedVector(residual);
00408
00409
00410 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00411 {
00412 DistributedVector::Stripe solution_stripe(solution_distributed, index_of_unknown);
00413 DistributedVector::Stripe residual_stripe(residual_distributed, index_of_unknown);
00414
00415 unsigned node_index = this->mDirichIterator->first->GetIndex();
00416
00417 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
00418
00419 if (solution_distributed.IsGlobalIndexLocal(node_index))
00420 {
00421 residual_stripe[node_index]=solution_stripe[node_index] - value;
00422 }
00423 this->mDirichIterator++;
00424 }
00425 solution_distributed.Restore();
00426 residual_distributed.Restore();
00427 }
00428 }
00429
00430 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00431 void BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyDirichletToNonlinearJacobian(Mat jacobian)
00432 {
00433 unsigned num_boundary_conditions = 0;
00434 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00435 {
00436 num_boundary_conditions += this->mpDirichletMap[index_of_unknown]->size();
00437 }
00438
00439 std::vector<unsigned> rows_to_zero(num_boundary_conditions);
00440
00441 unsigned counter=0;
00442 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00443 {
00444 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00445
00446 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00447 {
00448 unsigned node_index = this->mDirichIterator->first->GetIndex();
00449 rows_to_zero[counter++] = PROBLEM_DIM*node_index + index_of_unknown;
00450 this->mDirichIterator++;
00451 }
00452 }
00453 PetscMatTools::Finalise(jacobian);
00454 PetscMatTools::ZeroRowsWithValueOnDiagonal(jacobian, rows_to_zero, 1.0);
00455 }
00456
00457 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00458 bool BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::Validate(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh)
00459 {
00460 bool valid = true;
00461
00462 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00463 {
00464
00465 typename AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>::BoundaryElementIterator elt_iter
00466 = pMesh->GetBoundaryElementIteratorBegin();
00467 while (valid && elt_iter != pMesh->GetBoundaryElementIteratorEnd())
00468 {
00469 if (!HasNeumannBoundaryCondition(*elt_iter, index_of_unknown))
00470 {
00471
00472 for (unsigned i=0; i<(*elt_iter)->GetNumNodes(); i++)
00473 {
00474 if (!HasDirichletBoundaryCondition((*elt_iter)->GetNode(i)))
00475 {
00476 valid = false;
00477 }
00478 }
00479 }
00480 elt_iter++;
00481 }
00482 }
00483 return valid;
00484 }
00485
00486 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00487 double BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::GetNeumannBCValue(const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>* pSurfaceElement,
00488 const ChastePoint<SPACE_DIM>& rX,
00489 unsigned indexOfUnknown)
00490 {
00491 assert(indexOfUnknown < PROBLEM_DIM);
00492
00493
00494 if (mLastNeumannCondition[indexOfUnknown] == mpNeumannMap[indexOfUnknown]->end() ||
00495 mLastNeumannCondition[indexOfUnknown]->first != pSurfaceElement)
00496 {
00497 mLastNeumannCondition[indexOfUnknown] = mpNeumannMap[indexOfUnknown]->find(pSurfaceElement);
00498 }
00499 if (mLastNeumannCondition[indexOfUnknown] == mpNeumannMap[indexOfUnknown]->end())
00500 {
00501
00502 return 0.0;
00503 }
00504 else
00505 {
00506 return mLastNeumannCondition[indexOfUnknown]->second->GetValue(rX);
00507 }
00508 }
00509
00510 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00511 bool BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::HasNeumannBoundaryCondition(const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>* pSurfaceElement,
00512 unsigned indexOfUnknown)
00513 {
00514 assert(indexOfUnknown < PROBLEM_DIM);
00515
00516 mLastNeumannCondition[indexOfUnknown] = mpNeumannMap[indexOfUnknown]->find(pSurfaceElement);
00517
00518 return (mLastNeumannCondition[indexOfUnknown] != mpNeumannMap[indexOfUnknown]->end());
00519 }
00520
00521 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00522 bool BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::AnyNonZeroNeumannConditions()
00523 {
00524 bool ret = false;
00525 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00526 {
00527 if (mAnyNonZeroNeumannConditionsForUnknown[index_of_unknown] == true)
00528 {
00529 ret = true;
00530 }
00531 }
00532 return ret;
00533 }
00534
00535 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00536 typename BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::NeumannMapIterator BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::BeginNeumann()
00537 {
00538
00539 return mpNeumannMap[0]->begin();
00540 }
00541
00542 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00543 typename BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::NeumannMapIterator BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::EndNeumann()
00544 {
00545
00546 return mpNeumannMap[0]->end();
00547 }
00548
00549 #endif // _BOUNDARYCONDITIONSCONTAINERIMPLEMENTATION_HPP_