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