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 "BoundaryConditionsContainer.hpp"
00033 #include "DistributedVector.hpp"
00034 #include "ConstBoundaryCondition.hpp"
00035
00036 #include "PetscTools.hpp"
00037
00038
00039 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00040 BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::BoundaryConditionsContainer()
00041 : AbstractBoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>()
00042 {
00043 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00044 {
00045 mpNeumannMap[index_of_unknown] = new std::map< const BoundaryElement<ELEM_DIM-1, SPACE_DIM> *, const AbstractBoundaryCondition<SPACE_DIM>*>;
00046
00047 mAnyNonZeroNeumannConditionsForUnknown[index_of_unknown] = false;
00048 mLastNeumannCondition[index_of_unknown] = mpNeumannMap[index_of_unknown]->begin();
00049 }
00050
00051
00052 mpZeroBoundaryCondition = new ConstBoundaryCondition<SPACE_DIM>(0.0);
00053 }
00054
00055 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00056 BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::~BoundaryConditionsContainer()
00057 {
00058
00059
00060 std::set<const AbstractBoundaryCondition<SPACE_DIM>*> deleted_conditions;
00061 for (unsigned i=0; i<PROBLEM_DIM; i++)
00062 {
00063 NeumannMapIterator neumann_iterator = mpNeumannMap[i]->begin();
00064 while (neumann_iterator != mpNeumannMap[i]->end() )
00065 {
00066
00067 if (deleted_conditions.count(neumann_iterator->second) == 0)
00068 {
00069 deleted_conditions.insert(neumann_iterator->second);
00070
00071 if (neumann_iterator->second != mpZeroBoundaryCondition)
00072 {
00073 delete neumann_iterator->second;
00074 }
00075 }
00076 neumann_iterator++;
00077 }
00078 delete(mpNeumannMap[i]);
00079 }
00080
00081 delete mpZeroBoundaryCondition;
00082
00083 this->DeleteDirichletBoundaryConditions(deleted_conditions);
00084
00085 }
00086
00087 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00088 void BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::AddDirichletBoundaryCondition( const Node<SPACE_DIM> * pBoundaryNode,
00089 const AbstractBoundaryCondition<SPACE_DIM> * pBoundaryCondition,
00090 unsigned indexOfUnknown,
00091 bool checkIfBoundaryNode)
00092 {
00093 assert(indexOfUnknown < PROBLEM_DIM);
00094 if (checkIfBoundaryNode)
00095 {
00096 assert( pBoundaryNode->IsBoundaryNode());
00097 }
00098
00099 (*(this->mpDirichletMap[indexOfUnknown]))[pBoundaryNode] = pBoundaryCondition;
00100 }
00101
00102 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00103 void BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::AddNeumannBoundaryCondition( const BoundaryElement<ELEM_DIM-1, SPACE_DIM> * pBoundaryElement,
00104 const AbstractBoundaryCondition<SPACE_DIM> * pBoundaryCondition,
00105 unsigned indexOfUnknown)
00106 {
00107 assert(indexOfUnknown < PROBLEM_DIM);
00108
00109
00110 mAnyNonZeroNeumannConditionsForUnknown[indexOfUnknown] = true;
00111
00112 for (unsigned unknown=0; unknown<PROBLEM_DIM; unknown++)
00113 {
00114 if (unknown==indexOfUnknown)
00115 {
00116 (*(mpNeumannMap[indexOfUnknown]))[pBoundaryElement] = pBoundaryCondition;
00117 }
00118 else
00119 {
00120
00121 if ( mpNeumannMap[unknown]->find(pBoundaryElement)==mpNeumannMap[unknown]->end() )
00122 {
00123
00124 (*(mpNeumannMap[unknown]))[pBoundaryElement] = mpZeroBoundaryCondition;
00125 }
00126 }
00127 }
00128 }
00129
00130 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00131 void BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::DefineZeroDirichletOnMeshBoundary(AbstractTetrahedralMesh<ELEM_DIM,SPACE_DIM>* pMesh,
00132 unsigned indexOfUnknown)
00133 {
00134 this->DefineConstantDirichletOnMeshBoundary(pMesh, 0.0, indexOfUnknown);
00135 }
00136
00137 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00138 void BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::DefineConstantDirichletOnMeshBoundary(AbstractTetrahedralMesh<ELEM_DIM,SPACE_DIM>* pMesh,
00139 double value,
00140 unsigned indexOfUnknown)
00141 {
00142 assert(indexOfUnknown < PROBLEM_DIM);
00143
00144 assert(pMesh->GetNumBoundaryNodes() > 0);
00145
00146 ConstBoundaryCondition<SPACE_DIM> *p_boundary_condition =
00147 new ConstBoundaryCondition<SPACE_DIM>( value );
00148
00149 typename AbstractTetrahedralMesh<ELEM_DIM, SPACE_DIM>::BoundaryNodeIterator iter;
00150 iter = pMesh->GetBoundaryNodeIteratorBegin();
00151 while (iter != pMesh->GetBoundaryNodeIteratorEnd())
00152 {
00153 AddDirichletBoundaryCondition(*iter, p_boundary_condition, indexOfUnknown);
00154 iter++;
00155 }
00156 }
00157
00158 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00159 void BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::DefineZeroNeumannOnMeshBoundary(AbstractTetrahedralMesh<ELEM_DIM,SPACE_DIM> *pMesh,
00160 unsigned indexOfUnknown)
00161 {
00162 assert(indexOfUnknown < PROBLEM_DIM);
00163
00164 assert(pMesh->GetNumBoundaryElements() > 0);
00165 ConstBoundaryCondition<SPACE_DIM> *p_zero_boundary_condition =
00166 new ConstBoundaryCondition<SPACE_DIM>( 0.0 );
00167
00168 typename AbstractTetrahedralMesh<ELEM_DIM, SPACE_DIM>::BoundaryElementIterator iter;
00169 iter = pMesh->GetBoundaryElementIteratorBegin();
00170 while (iter != pMesh->GetBoundaryElementIteratorEnd())
00171 {
00172 AddNeumannBoundaryCondition(*iter, p_zero_boundary_condition, indexOfUnknown);
00173 iter++;
00174 }
00175
00176 mAnyNonZeroNeumannConditionsForUnknown[indexOfUnknown] = false;
00177 }
00208 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00209 void BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyDirichletToLinearProblem(LinearSystem& rLinearSystem,
00210 bool applyToMatrix)
00211 {
00212 if (applyToMatrix)
00213 {
00214
00215
00216
00217 VecDuplicate(rLinearSystem.rGetRhsVector(), &(rLinearSystem.rGetDirichletBoundaryConditionsVector()));
00218 VecZeroEntries(rLinearSystem.rGetDirichletBoundaryConditionsVector());
00219
00220
00221
00222 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00223 {
00224 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00225
00226 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00227 {
00228 unsigned node_index = this->mDirichIterator->first->GetIndex();
00229 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
00230
00231 unsigned row = PROBLEM_DIM*node_index + index_of_unknown;
00232 unsigned col = row;
00233
00234
00235
00236 Vec matrix_col;
00237 VecDuplicate(rLinearSystem.rGetRhsVector(), &matrix_col);
00238 VecZeroEntries(matrix_col);
00239
00240 rLinearSystem.AssembleFinalLinearSystem();
00241 Mat& r_mat = rLinearSystem.rGetLhsMatrix();
00242 MatGetColumnVector(r_mat, matrix_col, col);
00243
00244
00245 int indices[1] = {col};
00246 double zero[1] = {0.0};
00247 VecSetValues(matrix_col, 1, indices, zero, INSERT_VALUES);
00248
00249
00250
00251
00252
00253 VecAXPY(rLinearSystem.rGetDirichletBoundaryConditionsVector(), -value, matrix_col);
00254
00255 VecDestroy(matrix_col);
00256
00257
00258 rLinearSystem.ZeroMatrixRow(row);
00259 rLinearSystem.ZeroMatrixColumn(col);
00260 rLinearSystem.SetMatrixElement(row, row, 1);
00261
00262 this->mDirichIterator++;
00263
00264 }
00265 }
00266
00267 }
00268
00269
00270 if (rLinearSystem.rGetDirichletBoundaryConditionsVector())
00271 {
00272 VecAXPY(rLinearSystem.rGetRhsVector(), 1.0, rLinearSystem.rGetDirichletBoundaryConditionsVector());
00273 }
00274
00275
00276
00277 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00278 {
00279 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00280
00281 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00282 {
00283 unsigned node_index = this->mDirichIterator->first->GetIndex();
00284 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
00285
00286 unsigned row = PROBLEM_DIM*node_index + index_of_unknown;
00287
00288 rLinearSystem.SetRhsVectorElement(row, value);
00289
00290 this->mDirichIterator++;
00291 }
00292 }
00293 }
00294
00295 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00296 void BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyDirichletToNonlinearResidual(
00297 const Vec currentSolution,
00298 Vec residual,
00299 DistributedVectorFactory& rFactory)
00300 {
00301 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00302 {
00303 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00304
00305 DistributedVector solution_distributed = rFactory.CreateDistributedVector(currentSolution);
00306 DistributedVector residual_distributed = rFactory.CreateDistributedVector(residual);
00307
00308
00309 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00310 {
00311 DistributedVector::Stripe solution_stripe(solution_distributed, index_of_unknown);
00312 DistributedVector::Stripe residual_stripe(residual_distributed, index_of_unknown);
00313
00314 unsigned node_index = this->mDirichIterator->first->GetIndex();
00315
00316 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
00317
00318 if (solution_distributed.IsGlobalIndexLocal(node_index))
00319 {
00320 residual_stripe[node_index]=solution_stripe[node_index] - value;
00321 }
00322 this->mDirichIterator++;
00323 }
00324 solution_distributed.Restore();
00325 residual_distributed.Restore();
00326 }
00327 }
00328
00329 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00330 void BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::ApplyDirichletToNonlinearJacobian(Mat jacobian)
00331 {
00332 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00333 {
00334 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
00335 PetscInt irows, icols;
00336 double value;
00337 MatGetSize(jacobian, &irows, &icols);
00338 unsigned cols=icols;
00339
00340 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
00341 {
00342 unsigned node_index = this->mDirichIterator->first->GetIndex();
00343
00344 unsigned row_index = PROBLEM_DIM*node_index + index_of_unknown;
00345 assert(row_index<(unsigned)irows);
00346
00347 for (unsigned col_index=0; col_index<cols; col_index++)
00348 {
00349 value = (col_index == row_index) ? 1.0 : 0.0;
00350 MatSetValue(jacobian, row_index, col_index, value, INSERT_VALUES);
00351 }
00352 this->mDirichIterator++;
00353 }
00354 }
00355 }
00356
00357 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00358 bool BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::Validate(AbstractTetrahedralMesh<ELEM_DIM,SPACE_DIM>* pMesh)
00359 {
00360 bool valid = true;
00361
00362 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00363 {
00364
00365 typename AbstractTetrahedralMesh<ELEM_DIM,SPACE_DIM>::BoundaryElementIterator elt_iter
00366 = pMesh->GetBoundaryElementIteratorBegin();
00367 while (valid && elt_iter != pMesh->GetBoundaryElementIteratorEnd())
00368 {
00369 if (!HasNeumannBoundaryCondition(*elt_iter, index_of_unknown))
00370 {
00371
00372 for (unsigned i=0; i<(*elt_iter)->GetNumNodes(); i++)
00373 {
00374 if (!HasDirichletBoundaryCondition((*elt_iter)->GetNode(i)))
00375 {
00376 valid = false;
00377 }
00378 }
00379 }
00380 elt_iter++;
00381 }
00382 }
00383 return valid;
00384 }
00385
00386 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00387 double BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::GetNeumannBCValue(const BoundaryElement<ELEM_DIM-1,SPACE_DIM>* pSurfaceElement,
00388 const ChastePoint<SPACE_DIM>& rX,
00389 unsigned indexOfUnknown)
00390 {
00391 assert(indexOfUnknown < PROBLEM_DIM);
00392
00393
00394 if (mLastNeumannCondition[indexOfUnknown] == mpNeumannMap[indexOfUnknown]->end() ||
00395 mLastNeumannCondition[indexOfUnknown]->first != pSurfaceElement)
00396 {
00397 mLastNeumannCondition[indexOfUnknown] = mpNeumannMap[indexOfUnknown]->find(pSurfaceElement);
00398 }
00399 if (mLastNeumannCondition[indexOfUnknown] == mpNeumannMap[indexOfUnknown]->end())
00400 {
00401
00402 return 0.0;
00403 }
00404 else
00405 {
00406 return mLastNeumannCondition[indexOfUnknown]->second->GetValue(rX);
00407 }
00408 }
00409
00410 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00411 bool BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::HasNeumannBoundaryCondition(const BoundaryElement<ELEM_DIM-1,SPACE_DIM>* pSurfaceElement,
00412 unsigned indexOfUnknown)
00413 {
00414 assert(indexOfUnknown < PROBLEM_DIM);
00415
00416 mLastNeumannCondition[indexOfUnknown] = mpNeumannMap[indexOfUnknown]->find(pSurfaceElement);
00417
00418 return (mLastNeumannCondition[indexOfUnknown] != mpNeumannMap[indexOfUnknown]->end());
00419 }
00420
00421 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00422 bool BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::AnyNonZeroNeumannConditions()
00423 {
00424 bool ret = false;
00425 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00426 {
00427 if (mAnyNonZeroNeumannConditionsForUnknown[index_of_unknown] == true)
00428 {
00429 ret = true;
00430 }
00431 }
00432 return ret;
00433 }
00434
00435 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00436 typename BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::NeumannMapIterator BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::BeginNeumann()
00437 {
00438
00439 return mpNeumannMap[0]->begin();
00440 }
00441
00442 template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00443 typename BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::NeumannMapIterator BoundaryConditionsContainer<ELEM_DIM,SPACE_DIM,PROBLEM_DIM>::EndNeumann()
00444 {
00445
00446 return mpNeumannMap[0]->end();
00447 }
00448
00449 #endif // _BOUNDARYCONDITIONSCONTAINERIMPLEMENTATION_HPP_