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