Chaste Commit::675f9facbe008c5eacb9006feaeb6423206579ea
BoundaryConditionsContainerImplementation.hpp
1/*
2
3Copyright (c) 2005-2025, University of Oxford.
4All rights reserved.
5
6University of Oxford means the Chancellor, Masters and Scholars of the
7University of Oxford, having an administrative office at Wellington
8Square, Oxford OX1 2JD, UK.
9
10This file is part of Chaste.
11
12Redistribution and use in source and binary forms, with or without
13modification, are permitted provided that the following conditions are met:
14 * Redistributions of source code must retain the above copyright notice,
15 this list of conditions and the following disclaimer.
16 * Redistributions in binary form must reproduce the above copyright notice,
17 this list of conditions and the following disclaimer in the documentation
18 and/or other materials provided with the distribution.
19 * Neither the name of the University of Oxford nor the names of its
20 contributors may be used to endorse or promote products derived from this
21 software without specific prior written permission.
22
23THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33
34*/
35
36#ifndef _BOUNDARYCONDITIONSCONTAINERIMPLEMENTATION_HPP_
37#define _BOUNDARYCONDITIONSCONTAINERIMPLEMENTATION_HPP_
38
39#include "BoundaryConditionsContainer.hpp"
40#include "ConstBoundaryCondition.hpp"
41#include "DistributedVector.hpp"
42#include "Exception.hpp"
43#include "HeartEventHandler.hpp"
44#include "PetscMatTools.hpp"
45#include "PetscVecTools.hpp"
46#include "ReplicatableVector.hpp"
47
48template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
50 : AbstractBoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>(deleteConditions)
51{
52 mLoadedFromArchive = false;
53
54 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
55 {
56 mpNeumannMap[index_of_unknown] = new std::map< const BoundaryElement<ELEMENT_DIM-1, SPACE_DIM> *, const AbstractBoundaryCondition<SPACE_DIM>*>;
57
58 mAnyNonZeroNeumannConditionsForUnknown[index_of_unknown] = false;
59 mLastNeumannCondition[index_of_unknown] = mpNeumannMap[index_of_unknown]->begin();
60
61 mpPeriodicBcMap[index_of_unknown] = new std::map< const Node<SPACE_DIM> *, const Node<SPACE_DIM> * >;
62 }
63
64 // This zero boundary condition is only used in AddNeumannBoundaryCondition
66}
67
68template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
70{
71 // Keep track of what boundary condition objects we've deleted
72 std::set<const AbstractBoundaryCondition<SPACE_DIM>*> deleted_conditions;
73 for (unsigned i=0; i<PROBLEM_DIM; i++)
74 {
75 NeumannMapIterator neumann_iterator = mpNeumannMap[i]->begin();
76 while (neumann_iterator != mpNeumannMap[i]->end() )
77 {
78
79 if (deleted_conditions.count(neumann_iterator->second) == 0)
80 {
81 deleted_conditions.insert(neumann_iterator->second);
82 //Leave the zero boundary condition until last
83 if (neumann_iterator->second != mpZeroBoundaryCondition)
84 {
85 if (this->mDeleteConditions)
86 {
87 delete neumann_iterator->second;
88 }
89 }
90 }
91 neumann_iterator++;
92 }
93 delete(mpNeumannMap[i]);
94 delete(mpPeriodicBcMap[i]);
95 }
96
97 delete mpZeroBoundaryCondition;
98
99 if (this->mDeleteConditions)
100 {
101 this->DeleteDirichletBoundaryConditions(deleted_conditions);
102 }
103}
104
105template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
107 const AbstractBoundaryCondition<SPACE_DIM>* pBoundaryCondition,
108 unsigned indexOfUnknown,
109 bool checkIfBoundaryNode)
111 assert(indexOfUnknown < PROBLEM_DIM);
112 if (checkIfBoundaryNode)
113 {
114 assert(pBoundaryNode->IsBoundaryNode());
115 }
117 (*(this->mpDirichletMap[indexOfUnknown]))[pBoundaryNode] = pBoundaryCondition;
118}
119
120template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
122 const Node<SPACE_DIM>* pNode2)
123{
124 assert(pNode1->IsBoundaryNode());
125 assert(pNode2->IsBoundaryNode());
126
127 // will assume the periodic BC is to be applied to ALL unknowns, can't really imagine a
128 // situation where this isn't going to be true. If necessary can easily change this method
129 // to take in the index of the unknown
130 for (unsigned i=0; i<PROBLEM_DIM; i++)
131 {
132 (*(this->mpPeriodicBcMap[i]))[pNode1] = pNode2;
133 }
134}
135
136template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
138 const AbstractBoundaryCondition<SPACE_DIM> * pBoundaryCondition,
139 unsigned indexOfUnknown)
140{
141 assert(indexOfUnknown < PROBLEM_DIM);
142
143 /*
144 * If this condition is constant, we can test whether it is zero.
145 * Otherwise we assume that this could be a non-zero boundary condition.
146 */
147 const ConstBoundaryCondition<SPACE_DIM>* p_const_cond = dynamic_cast<const ConstBoundaryCondition<SPACE_DIM>*>(pBoundaryCondition);
148 if (p_const_cond)
149 {
150 if (p_const_cond->GetValue(pBoundaryElement->GetNode(0)->GetPoint()) != 0.0)
151 {
152 mAnyNonZeroNeumannConditionsForUnknown[indexOfUnknown] = true;
154 }
155 else
156 {
157 mAnyNonZeroNeumannConditionsForUnknown[indexOfUnknown] = true;
158 }
159
160 for (unsigned unknown=0; unknown<PROBLEM_DIM; unknown++)
161 {
162 if (unknown == indexOfUnknown)
163 {
164 (*(mpNeumannMap[indexOfUnknown]))[pBoundaryElement] = pBoundaryCondition;
166 else
167 {
168 // If can't find pBoundaryElement in map[unknown]
169 if (mpNeumannMap[unknown]->find(pBoundaryElement)==mpNeumannMap[unknown]->end())
170 {
171 // Add zero bc to other unknowns (so all maps are in sync)
172 (*(mpNeumannMap[unknown]))[pBoundaryElement] = mpZeroBoundaryCondition;
173 }
174 }
175 }
177
178template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
180 unsigned indexOfUnknown)
181{
182 this->DefineConstantDirichletOnMeshBoundary(pMesh, 0.0, indexOfUnknown);
183}
184
185template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
187 double value,
188 unsigned indexOfUnknown)
189{
190 assert(indexOfUnknown < PROBLEM_DIM);
191 // In applying a condition to the boundary, we need to be sure that the boundary exists
192 assert(PetscTools::ReplicateBool( pMesh->GetNumBoundaryNodes() > 0 ) );
193
195
197 iter = pMesh->GetBoundaryNodeIteratorBegin();
198 while (iter != pMesh->GetBoundaryNodeIteratorEnd())
199 {
200 AddDirichletBoundaryCondition(*iter, p_boundary_condition, indexOfUnknown);
201 iter++;
202 }
203}
204
205template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
207 unsigned indexOfUnknown)
208{
209 assert(indexOfUnknown < PROBLEM_DIM);
210
211 // In applying a condition to the boundary, we need to be sure that the boundary exists
212 assert(pMesh->GetNumBoundaryElements() > 0);
213 ConstBoundaryCondition<SPACE_DIM>* p_zero_boundary_condition = new ConstBoundaryCondition<SPACE_DIM>( 0.0 );
214
216 iter = pMesh->GetBoundaryElementIteratorBegin();
217 while (iter != pMesh->GetBoundaryElementIteratorEnd())
218 {
219 AddNeumannBoundaryCondition(*iter, p_zero_boundary_condition, indexOfUnknown);
220 iter++;
221 }
222
223 mAnyNonZeroNeumannConditionsForUnknown[indexOfUnknown] = false;
224}
225
256template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
258 LinearSystem& rLinearSystem,
259 bool applyToMatrix,
260 bool applyToRhsVector)
261{
262 HeartEventHandler::BeginEvent(HeartEventHandler::DIRICHLET_BCS);
263
264 if (applyToMatrix)
265 {
266 if (!this->HasDirichletBoundaryConditions())
268 // Short-circuit the replication if there are no conditions
269 HeartEventHandler::EndEvent(HeartEventHandler::DIRICHLET_BCS);
270 return;
271 }
272
273 bool matrix_is_symmetric = rLinearSystem.IsMatrixSymmetric();
274
275 if (matrix_is_symmetric)
276 {
277 /*
278 * Modifications to the RHS are stored in the Dirichlet boundary
279 * conditions vector. This is done so that they can be reapplied
280 * at each time step.
281 * Make a new vector to store the Dirichlet offsets in.
282 */
283 Vec& r_bcs_vec = rLinearSystem.rGetDirichletBoundaryConditionsVector();
284 if (!r_bcs_vec)
285 {
286 VecDuplicate(rLinearSystem.rGetRhsVector(), &r_bcs_vec);
287 }
288 PetscVecTools::Zero(r_bcs_vec);
289 /*
290 * If the matrix is symmetric, calls to GetMatrixRowDistributed()
291 * require the matrix to be in assembled state. Otherwise we can
292 * defer it.
293 */
294 rLinearSystem.AssembleFinalLinearSystem();
295 }
296
297 // Work out where we're setting Dirichlet boundary conditions *everywhere*, not just those locally known
298 ReplicatableVector dirichlet_conditions(rLinearSystem.GetSize());
299 unsigned lo, hi;
300 {
301 PetscInt lo_s, hi_s;
302 rLinearSystem.GetOwnershipRange(lo_s, hi_s);
303 lo = lo_s; hi = hi_s;
304 }
305 // Initialise all local entries to DBL_MAX, i.e. don't know if there's a condition
306 for (unsigned i=lo; i<hi; i++)
307 {
308 dirichlet_conditions[i] = DBL_MAX;
310 // Now fill in the ones we know
311 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
312 {
313 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
315 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
316 {
317 unsigned node_index = this->mDirichIterator->first->GetIndex();
318 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
319 assert(value != DBL_MAX);
320
321 unsigned row = PROBLEM_DIM*node_index + index_of_unknown;
322 dirichlet_conditions[row] = value;
323
324 this->mDirichIterator++;
325 }
326 }
327
328 // And replicate
329 dirichlet_conditions.Replicate(lo, hi);
330
331 // Which rows have conditions?
332 std::vector<unsigned> rows_to_zero;
333 for (unsigned i=0; i<dirichlet_conditions.GetSize(); i++)
334 {
335 if (dirichlet_conditions[i] != DBL_MAX)
336 {
337 rows_to_zero.push_back(i);
338 }
339 }
340
341 if (matrix_is_symmetric)
342 {
343 // Modify the matrix columns
344 for (unsigned i=0; i<rows_to_zero.size(); i++)
345 {
346 unsigned col = rows_to_zero[i];
347 double minus_value = -dirichlet_conditions[col];
348
349 /*
350 * Get a vector which will store the column of the matrix (column d,
351 * where d is the index of the row (and column) to be altered for the
352 * boundary condition. Since the matrix is symmetric when get row
353 * number "col" and treat it as a column. PETSc uses compressed row
354 * format and therefore getting rows is far more efficient than getting
355 * columns.
356 */
357 Vec matrix_col = rLinearSystem.GetMatrixRowDistributed(col);
358
359 // Zero the correct entry of the column
360 PetscVecTools::SetElement(matrix_col, col, 0.0);
361
362 /*
363 * Set up the RHS Dirichlet boundary conditions vector.
364 * Assuming one boundary at the zeroth node (x_0 = value), this is equal to
365 * -value*[0 a_21 a_31 .. a_N1]
366 * and will be added to the RHS.
367 */
368 PetscVecTools::Finalise(matrix_col); //PETSc 3.19 doesn't like mixing "set" mode (above) with "add" mode
369 PetscVecTools::AddScaledVector(rLinearSystem.rGetDirichletBoundaryConditionsVector(), matrix_col, minus_value);
370 PetscTools::Destroy(matrix_col);
371 }
372 }
373
374 /*
375 * Now zero the appropriate rows and columns of the matrix. If the matrix
376 * is symmetric we apply the boundary conditions in a way the symmetry isn't
377 * lost (rows and columns). If not only the row is zeroed.
378 */
379 if (matrix_is_symmetric)
380 {
381 rLinearSystem.ZeroMatrixRowsAndColumnsWithValueOnDiagonal(rows_to_zero, 1.0);
382 }
383 else
384 {
385 rLinearSystem.ZeroMatrixRowsWithValueOnDiagonal(rows_to_zero, 1.0);
386 }
387 }
388
389 if (applyToRhsVector)
390 {
391 // Apply the RHS boundary conditions modification if required.
392 if (rLinearSystem.rGetDirichletBoundaryConditionsVector())
393 {
395 }
396
397 /*
398 * Apply the actual boundary condition to the RHS, note this must be
399 * done after the modification to the RHS vector.
400 */
401 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
402 {
403 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
404
405 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
406 {
407 unsigned node_index = this->mDirichIterator->first->GetIndex();
408 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
409
410 unsigned row = PROBLEM_DIM*node_index + index_of_unknown;
411
412 rLinearSystem.SetRhsVectorElement(row, value);
413
414 this->mDirichIterator++;
415 }
416 }
417 }
418
419 HeartEventHandler::EndEvent(HeartEventHandler::DIRICHLET_BCS);
420}
421
422template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
424 bool applyToMatrix,
425 bool applyToRhsVector)
426{
427 bool has_periodic_bcs = false;
428 for (unsigned i=0; i<PROBLEM_DIM; i++)
429 {
430 if (!mpPeriodicBcMap[i]->empty())
431 {
432 has_periodic_bcs = true;
433 break;
434 }
435 }
436
437 EXCEPT_IF_NOT(has_periodic_bcs);
438
439 if (applyToMatrix)
440 {
441 std::vector<unsigned> rows_to_zero;
442 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
443 {
444 for (typename std::map< const Node<SPACE_DIM> *, const Node<SPACE_DIM> * >::const_iterator iter = mpPeriodicBcMap[index_of_unknown]->begin();
445 iter != mpPeriodicBcMap[index_of_unknown]->end();
446 ++iter)
447 {
448 unsigned node_index_1 = iter->first->GetIndex();
449 unsigned row_index_1 = PROBLEM_DIM*node_index_1 + index_of_unknown;
450 rows_to_zero.push_back(row_index_1);
451 }
452 }
453
454 rLinearSystem.ZeroMatrixRowsWithValueOnDiagonal(rows_to_zero, 1.0);
455
456 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
457 {
458 for (typename std::map< const Node<SPACE_DIM> *, const Node<SPACE_DIM> * >::const_iterator iter = mpPeriodicBcMap[index_of_unknown]->begin();
459 iter != mpPeriodicBcMap[index_of_unknown]->end();
460 ++iter)
461 {
462 unsigned node_index_1 = iter->first->GetIndex();
463 unsigned node_index_2 = iter->second->GetIndex();
464
465 unsigned mat_index1 = PROBLEM_DIM*node_index_1 + index_of_unknown;
466 unsigned mat_index2 = PROBLEM_DIM*node_index_2 + index_of_unknown;
467 PetscMatTools::SetElement(rLinearSystem.rGetLhsMatrix(), mat_index1, mat_index2, -1.0);
468 }
469 }
470 }
471
472 if (applyToRhsVector)
473 {
474 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
475 {
476 for (typename std::map< const Node<SPACE_DIM> *, const Node<SPACE_DIM> * >::const_iterator iter = mpPeriodicBcMap[index_of_unknown]->begin();
477 iter != mpPeriodicBcMap[index_of_unknown]->end();
478 ++iter)
479 {
480 unsigned node_index = iter->first->GetIndex();
481 unsigned row_index = PROBLEM_DIM*node_index + index_of_unknown;
482 rLinearSystem.SetRhsVectorElement(row_index, 0.0);
483 }
484 }
485 }
486}
487
488template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
490 const Vec currentSolution,
491 Vec residual,
492 DistributedVectorFactory& rFactory)
493{
494 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
495 {
496 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
497
498 DistributedVector solution_distributed = rFactory.CreateDistributedVector(currentSolution, true /*Read-only*/);
499 DistributedVector residual_distributed = rFactory.CreateDistributedVector(residual);
500
501
502 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
503 {
504 DistributedVector::Stripe solution_stripe(solution_distributed, index_of_unknown);
505 DistributedVector::Stripe residual_stripe(residual_distributed, index_of_unknown);
506
507 unsigned node_index = this->mDirichIterator->first->GetIndex();
508
509 double value = this->mDirichIterator->second->GetValue(this->mDirichIterator->first->GetPoint());
510
511 if (solution_distributed.IsGlobalIndexLocal(node_index))
512 {
513 residual_stripe[node_index]=solution_stripe[node_index] - value;
514 }
515 this->mDirichIterator++;
516 }
517 // Don't restore the read-only one: solution_distributed.Restore();
518 residual_distributed.Restore();
519 }
520}
521
522template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
524{
525 unsigned num_boundary_conditions = 0;
526 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
527 {
528 num_boundary_conditions += this->mpDirichletMap[index_of_unknown]->size();
529 }
530
531 std::vector<unsigned> rows_to_zero(num_boundary_conditions);
532
533 unsigned counter=0;
534 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
535 {
536 this->mDirichIterator = this->mpDirichletMap[index_of_unknown]->begin();
537
538 while (this->mDirichIterator != this->mpDirichletMap[index_of_unknown]->end() )
539 {
540 unsigned node_index = this->mDirichIterator->first->GetIndex();
541 rows_to_zero[counter++] = PROBLEM_DIM*node_index + index_of_unknown;
542 this->mDirichIterator++;
543 }
544 }
545 PetscMatTools::Finalise(jacobian);
546 PetscMatTools::ZeroRowsWithValueOnDiagonal(jacobian, rows_to_zero, 1.0);
547}
548
549template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
551{
552 bool valid = true;
553
554 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
555 {
556 // Iterate over surface elements
559 while (valid && elt_iter != pMesh->GetBoundaryElementIteratorEnd())
560 {
561 if (!HasNeumannBoundaryCondition(*elt_iter, index_of_unknown))
562 {
563 // Check for Dirichlet conditions on this element's nodes
564 for (unsigned i=0; i<(*elt_iter)->GetNumNodes(); i++)
565 {
566 if (!this->HasDirichletBoundaryCondition((*elt_iter)->GetNode(i)))
567 {
568 valid = false;
569 }
570 }
571 }
572 elt_iter++;
573 }
574 }
575 return valid;
576}
577
578template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
580 const ChastePoint<SPACE_DIM>& rX,
581 unsigned indexOfUnknown)
582{
583 assert(indexOfUnknown < PROBLEM_DIM);
584
585 // Did we see this condition on the last search we did?
586 if (mLastNeumannCondition[indexOfUnknown] == mpNeumannMap[indexOfUnknown]->end() ||
587 mLastNeumannCondition[indexOfUnknown]->first != pSurfaceElement)
588 {
589 mLastNeumannCondition[indexOfUnknown] = mpNeumannMap[indexOfUnknown]->find(pSurfaceElement);
590 }
591 if (mLastNeumannCondition[indexOfUnknown] == mpNeumannMap[indexOfUnknown]->end())
592 {
593 // No Neumann condition is equivalent to a zero Neumann condition
594 return 0.0;
595 }
596 else
597 {
598 return mLastNeumannCondition[indexOfUnknown]->second->GetValue(rX);
599 }
600}
601
602template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
604 unsigned indexOfUnknown)
605{
606 assert(indexOfUnknown < PROBLEM_DIM);
607
608 mLastNeumannCondition[indexOfUnknown] = mpNeumannMap[indexOfUnknown]->find(pSurfaceElement);
609
610 return (mLastNeumannCondition[indexOfUnknown] != mpNeumannMap[indexOfUnknown]->end());
611}
612
613template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
615{
616 bool ret = false;
617 for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
618 {
619 if (mAnyNonZeroNeumannConditionsForUnknown[index_of_unknown] == true)
620 {
621 ret = true;
622 }
623 }
624 return ret;
625}
626
627template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
629{
630 // [0] is ok as all maps will be in sync due to the way ApplyNeumannBoundaryCondition works
631 return mpNeumannMap[0]->begin();
632}
633
634template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
636{
637 // [0] is ok as all maps will be in sync due to the way ApplyNeumannBoundaryCondition works
638 return mpNeumannMap[0]->end();
639}
640
641#endif // _BOUNDARYCONDITIONSCONTAINERIMPLEMENTATION_HPP_
#define EXCEPT_IF_NOT(test)
Node< SPACE_DIM > * GetNode(unsigned localIndex) const
BoundaryNodeIterator GetBoundaryNodeIteratorBegin() const
std::vector< Node< SPACE_DIM > * >::const_iterator BoundaryNodeIterator
BoundaryNodeIterator GetBoundaryNodeIteratorEnd() const
unsigned GetNumBoundaryNodes() const
virtual unsigned GetNumBoundaryElements() const
std::vector< BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > * >::const_iterator BoundaryElementIterator
BoundaryElementIterator GetBoundaryElementIteratorBegin() const
BoundaryElementIterator GetBoundaryElementIteratorEnd() const
std::map< const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > *, const AbstractBoundaryCondition< SPACE_DIM > * > * mpNeumannMap[PROBLEM_DIM]
NeumannMapIterator mLastNeumannCondition[PROBLEM_DIM]
void AddDirichletBoundaryCondition(const Node< SPACE_DIM > *pBoundaryNode, const AbstractBoundaryCondition< SPACE_DIM > *pBoundaryCondition, unsigned indexOfUnknown=0, bool checkIfBoundaryNode=true)
double GetNeumannBCValue(const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > *pSurfaceElement, const ChastePoint< SPACE_DIM > &rX, unsigned indexOfUnknown=0)
void ApplyDirichletToLinearProblem(LinearSystem &rLinearSystem, bool applyToMatrix=true, bool applyToRhsVector=true)
void AddPeriodicBoundaryCondition(const Node< SPACE_DIM > *pNode1, const Node< SPACE_DIM > *pNode2)
bool HasNeumannBoundaryCondition(const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > *pSurfaceElement, unsigned indexOfUnknown=0)
bool Validate(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
void DefineZeroDirichletOnMeshBoundary(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, unsigned indexOfUnknown=0)
std::map< const Node< SPACE_DIM > *, const Node< SPACE_DIM > * > * mpPeriodicBcMap[PROBLEM_DIM]
void ApplyDirichletToNonlinearResidual(const Vec currentSolution, Vec residual, DistributedVectorFactory &rFactory)
void DefineZeroNeumannOnMeshBoundary(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, unsigned indexOfUnknown=0)
std::map< constBoundaryElement< ELEMENT_DIM-1, SPACE_DIM > *, constAbstractBoundaryCondition< SPACE_DIM > * >::const_iterator NeumannMapIterator
ConstBoundaryCondition< SPACE_DIM > * mpZeroBoundaryCondition
void DefineConstantDirichletOnMeshBoundary(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, double value, unsigned indexOfUnknown=0)
void AddNeumannBoundaryCondition(const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > *pBoundaryElement, const AbstractBoundaryCondition< SPACE_DIM > *pBoundaryCondition, unsigned indexOfUnknown=0)
void ApplyPeriodicBcsToLinearProblem(LinearSystem &rLinearSystem, bool applyToMatrix=true, bool applyToRhsVector=true)
double GetValue(const ChastePoint< SPACE_DIM > &rX) const
DistributedVector CreateDistributedVector(Vec vec, bool readOnly=false)
bool IsGlobalIndexLocal(unsigned globalIndex)
Vec & rGetDirichletBoundaryConditionsVector()
void AssembleFinalLinearSystem()
void GetOwnershipRange(PetscInt &lo, PetscInt &hi)
Mat & rGetLhsMatrix()
bool IsMatrixSymmetric()
void ZeroMatrixRowsWithValueOnDiagonal(std::vector< unsigned > &rRows, double diagonalValue)
unsigned GetSize() const
void ZeroMatrixRowsAndColumnsWithValueOnDiagonal(std::vector< unsigned > &rRowColIndices, double diagonalValue)
void SetRhsVectorElement(PetscInt row, double value)
Vec GetMatrixRowDistributed(unsigned rowIndex)
Vec & rGetRhsVector()
Definition Node.hpp:59
bool IsBoundaryNode() const
Definition Node.cpp:164
static void ZeroRowsWithValueOnDiagonal(Mat matrix, std::vector< unsigned > &rRows, double diagonalValue)
static void Finalise(Mat matrix)
static void SetElement(Mat matrix, PetscInt row, PetscInt col, double value)
static void Destroy(Vec &rVec)
static bool ReplicateBool(bool flag)
static void AddScaledVector(Vec y, Vec x, double scaleFactor)
static void Finalise(Vec vector)
static void SetElement(Vec vector, PetscInt row, double value)
static void Zero(Vec vector)