Chaste  Release::3.4
ElectrodesStimulusFactory.cpp
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 #include "ElectrodesStimulusFactory.hpp"
37 #include "DistributedTetrahedralMesh.hpp"
38 #include "IsNan.hpp"
39 #include "HeartConfig.hpp"
40 #include "GaussianQuadratureRule.hpp"
41 #include "RegularStimulus.hpp"
42 
43 template<unsigned DIM>
45  std::vector<double>& rStimulusMagnitudes,
46  std::vector<double>& rDurations,
47  std::vector<double>& rPeriods,
48  std::vector<double>& rStarts,
49  std::vector<double>& rEnds)
50  : mrElectrodePairs(rElectrodePairs),
51  mrMagnitudes(rStimulusMagnitudes),
52  mrDurations(rDurations),
53  mrPeriods(rPeriods),
54  mrStarts(rStarts),
55  mrEnds(rEnds),
56  mGroundSecondElectrode(false)
57 {
58 
59  if ( ( rElectrodePairs.size() != rStimulusMagnitudes.size() ) ||
60  ( rElectrodePairs.size() != rDurations.size() ) ||
61  ( rElectrodePairs.size() != rPeriods.size() ) ||
62  ( rElectrodePairs.size() != rStarts.size() ) ||
63  ( rElectrodePairs.size() != rEnds.size() ) )
64  {
65  EXCEPTION ("Vector of electrode pairs and vector of stimulation paremeters must have the same size");
66  }
67 
70 }
71 
72 template<unsigned DIM>
74 {
75 }
76 
77 template<unsigned DIM>
79 {
80  std::vector<unsigned> nodes_in_all_electrodes;
81  for (unsigned global_node_index = 0; global_node_index < this->mpMesh->GetNumNodes(); global_node_index++)
82  {
83  if (this->mpMesh->GetDistributedVectorFactory()->IsGlobalIndexLocal(global_node_index))
84  {
85  for (unsigned pair_index = 0; pair_index <mrElectrodePairs.size(); pair_index++)
86  {
87  if ( mrElectrodePairs[pair_index].first->DoesContain( this->mpMesh->GetNode(global_node_index)->GetPoint() ) )
88  {
89  nodes_in_all_electrodes.push_back( this->mpMesh->GetNode(global_node_index)->GetIndex() );
90  }
91  if ( mrElectrodePairs[pair_index].second->DoesContain( this->mpMesh->GetNode(global_node_index)->GetPoint() ) )
92  {
93  nodes_in_all_electrodes.push_back( this->mpMesh->GetNode(global_node_index)->GetIndex() );
94  }
95  }
96  }
97  }
99  for (unsigned node_index = 0; node_index < nodes_in_all_electrodes.size(); node_index++)
100  {
101  unsigned number_of_hits = 0;
102  for (unsigned node_to_check = 0; node_to_check < nodes_in_all_electrodes.size(); node_to_check++)
103  {
104  if (nodes_in_all_electrodes[node_index] == nodes_in_all_electrodes[node_to_check] )
105  {
106  number_of_hits++;
107  }
108  }
109  if (number_of_hits>1)
110  {
111  EXCEPTION("Two or more electrodes intersect with each other");
112  }
113  }
114 }
115 
116 template<unsigned DIM>
118 {
119  mGroundSecondElectrode = grounded;
120 }
121 
122 template<unsigned DIM>
124 {
125  assert(this->mpMesh!=NULL);
126  try
127  {
128  CheckForElectrodesIntersection();
129  }
130  catch (Exception &e)
131  {
132  PetscTools::ReplicateException(true); //Electrodes intersect
133  throw e;
134  }
136 
137  for (unsigned pair_index = 0; pair_index < mrElectrodePairs.size(); pair_index++)
138  {
139 
140  if (!mGroundSecondElectrode)//no grounding of second electrode
141  {
142  //compute the two fluxes for this pair
143  double flux_electrode_1 = ComputeElectrodeTotalFlux(mrElectrodePairs[pair_index].first, mMagnitudesElectrode1[pair_index]);
144  double flux_electrode_2 = ComputeElectrodeTotalFlux(mrElectrodePairs[pair_index].second, mMagnitudesElectrode2[pair_index]);
145 
146  //Scale the magnitude of the second electrode for this pair
147  mMagnitudesElectrode2[pair_index] = -mMagnitudesElectrode1[pair_index]*flux_electrode_1/flux_electrode_2;
148 
149  // Paranoia.
150  assert( flux_electrode_2 != 0.0);
151  assert( ! std::isnan(mMagnitudesElectrode2[pair_index]));
152 
153  }
154  else//second electrode is grounded
155  {
156  this->mGroundedRegions.push_back( mrElectrodePairs[pair_index].second );
157  }
158  }
159 }
160 
161 template<unsigned DIM>
162 boost::shared_ptr<AbstractStimulusFunction> ElectrodesStimulusFactory<DIM>::CreateStimulusForNode(Node<DIM>* pNode)
163 {
164  boost::shared_ptr<RegularStimulus> p_stimulus;
165  for (unsigned pair_index = 0; pair_index < mrElectrodePairs.size(); pair_index++)
166  {
167  if (mrElectrodePairs[pair_index].first->DoesContain(pNode->GetPoint()) )
168  {
169  p_stimulus.reset( new RegularStimulus(mMagnitudesElectrode1[pair_index], mrDurations[pair_index], mrPeriods[pair_index], mrStarts[pair_index], mrEnds[pair_index]));
170 
171  }
172  else if (mrElectrodePairs[pair_index].second->DoesContain(pNode->GetPoint()) )
173  {
174  p_stimulus.reset ( new RegularStimulus(mMagnitudesElectrode2[pair_index], mrDurations[pair_index], mrPeriods[pair_index], mrStarts[pair_index], mrEnds[pair_index]));
175 
176  }
177  else//no stimulus here
178  {
179  p_stimulus.reset ( new RegularStimulus(0.0, mrDurations[pair_index], mrPeriods[pair_index], mrStarts[pair_index], mrEnds[pair_index]));
180  }
181  }
182  return p_stimulus;
183 }
184 
185 template<unsigned DIM>
187 {
189 
190  //the basis functions
191  c_vector<double, DIM+1> phi;
192 
193  double total_electrode_flux = 0.0;
194  double ret;
195 
196  for (typename AbstractTetrahedralMesh<DIM,DIM>::NodeIterator node_iter=this->mpMesh->GetNodeIteratorBegin();
197  node_iter != this->mpMesh->GetNodeIteratorEnd();
198  ++node_iter)
199  {
200  if ( pRegion->DoesContain( (*node_iter).GetPoint() ) )
201  {
202  unsigned node_index = node_iter->GetIndex();
203  assert(node_index < this->mpMesh->GetNumNodes());
204  double contribution_of_this_node = 0.0;
205  //loop over the elements where this node is contained
206  for(std::set<unsigned>::iterator iter = this->mpMesh->GetNode(node_index)->rGetContainingElementIndices().begin();
207  iter != this->mpMesh->GetNode(node_index)->rGetContainingElementIndices().end();
208  ++iter)
209  {
210  Element<DIM,DIM>* p_element = this->mpMesh->GetElement(*iter);
211 
212  /*Determine jacobian for this element*/
213  c_matrix<double, DIM, DIM> jacobian;
214  c_matrix<double, DIM, DIM> inverse_jacobian;//unused here, but needed for the function below
215  double jacobian_determinant;
216  this->mpMesh->GetInverseJacobianForElement(p_element->GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
217 
218  double contribution_of_this_element = 0.0;//...to this node
219  // loop over Gauss points
220  for (unsigned quad_index=0; quad_index < pQuadRule->GetNumQuadPoints(); quad_index++)
221  {
222  const ChastePoint<DIM>& quad_point = pQuadRule->rGetQuadPoint(quad_index);
223  BasisFunction::ComputeBasisFunctions(quad_point, phi);
224 
225  double interpolated_stimulus = 0.0;
226  //loop over nodes in this element
227  for (unsigned node_index_in_element = 0; node_index_in_element < p_element->GetNumNodes(); node_index_in_element++)
228  {
229  //const Node<DIM>* p_node = p_element->GetNode(node_index_in_element);
230  assert(p_element->GetNumNodes() == DIM+1);
231  interpolated_stimulus += stimulusMagnitude*phi(node_index_in_element);
232  contribution_of_this_element += interpolated_stimulus * phi(node_index_in_element) * jacobian_determinant * pQuadRule->GetWeight(quad_index);
233  }
234 
235  }/*end of loop over gauss points*/
236  contribution_of_this_node += contribution_of_this_element;
237 
238  }/*end of loop over elements where the node is contained*/
239  total_electrode_flux += contribution_of_this_node;
240  }/* end of if that checks if node is in the electrode*/
241  }/* end of loop over nodes in the mesh*/
242 #ifndef NDEBUG
243  int mpi_ret = MPI_Allreduce(&total_electrode_flux, &ret, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
244  assert(mpi_ret == MPI_SUCCESS);
245 #else
246  MPI_Allreduce(&total_electrode_flux, &ret, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
247 #endif
248 
249  //clear up memory
250  delete pQuadRule;
251 
252  assert(ret < DBL_MAX);
253  return ret;
254 }
255 
257 // Explicit instantiation
259 
260 template class ElectrodesStimulusFactory<1>;
261 template class ElectrodesStimulusFactory<2>;
262 template class ElectrodesStimulusFactory<3>;
263 
std::vector< double > mMagnitudesElectrode2
Definition: Node.hpp:58
ElectrodesStimulusFactory(std::vector< std::pair< AbstractChasteRegion< DIM > *, AbstractChasteRegion< DIM > * > > &rElectrodePairs, std::vector< double > &rStimulusMagnitudes, std::vector< double > &rDurations, std::vector< double > &rPeriods, std::vector< double > &rStarts, std::vector< double > &rEnds)
static void Barrier(const std::string callerId="")
Definition: PetscTools.cpp:134
#define EXCEPTION(message)
Definition: Exception.hpp:143
std::vector< double > & mrMagnitudes
double ComputeElectrodeTotalFlux(AbstractChasteRegion< DIM > *pRegion, double stimulusMagnitude)
virtual bool DoesContain(const ChastePoint< SPACE_DIM > &rPointToCheck) const =0
unsigned GetNumQuadPoints() const
std::vector< double > mMagnitudesElectrode1
unsigned GetNumNodes() const
static void ReplicateException(bool flag)
Definition: PetscTools.cpp:198
NodeIterator GetNodeIteratorBegin(bool skipDeletedNodes=true)
unsigned GetIndex() const
ChastePoint< SPACE_DIM > GetPoint() const
Definition: Node.cpp:134
double GetWeight(unsigned index) const
boost::shared_ptr< AbstractStimulusFunction > CreateStimulusForNode(Node< DIM > *pNode)
const ChastePoint< ELEMENT_DIM > & rGetQuadPoint(unsigned index) const