36 #include "ElectrodesStimulusFactory.hpp" 37 #include "DistributedTetrahedralMesh.hpp" 39 #include "HeartConfig.hpp" 40 #include "GaussianQuadratureRule.hpp" 41 #include "RegularStimulus.hpp" 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),
56 mGroundSecondElectrode(false)
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()))
64 EXCEPTION (
"Vector of electrode pairs and vector of stimulation paremeters must have the same size");
71 template<
unsigned DIM>
76 template<
unsigned DIM>
79 std::vector<unsigned> nodes_in_all_electrodes;
80 for (
unsigned global_node_index = 0; global_node_index < this->
mpMesh->
GetNumNodes(); global_node_index++)
84 for (
unsigned pair_index = 0; pair_index <
mrElectrodePairs.size(); pair_index++)
86 if (
mrElectrodePairs[pair_index].first->DoesContain( this->mpMesh->GetNode(global_node_index)->GetPoint()))
88 nodes_in_all_electrodes.push_back( this->
mpMesh->
GetNode(global_node_index)->GetIndex() );
90 if (
mrElectrodePairs[pair_index].second->DoesContain( this->mpMesh->GetNode(global_node_index)->GetPoint()))
92 nodes_in_all_electrodes.push_back( this->
mpMesh->
GetNode(global_node_index)->GetIndex() );
98 for (
unsigned node_index = 0; node_index < nodes_in_all_electrodes.size(); node_index++)
100 unsigned number_of_hits = 0;
101 for (
unsigned node_to_check = 0; node_to_check < nodes_in_all_electrodes.size(); node_to_check++)
103 if (nodes_in_all_electrodes[node_index] == nodes_in_all_electrodes[node_to_check] )
108 if (number_of_hits>1)
110 EXCEPTION(
"Two or more electrodes intersect with each other");
115 template<
unsigned DIM>
121 template<
unsigned DIM>
124 assert(this->
mpMesh!=NULL);
136 for (
unsigned pair_index = 0; pair_index <
mrElectrodePairs.size(); pair_index++)
149 assert( flux_electrode_2 != 0.0);
160 template<
unsigned DIM>
163 boost::shared_ptr<RegularStimulus> p_stimulus;
164 for (
unsigned pair_index = 0; pair_index <
mrElectrodePairs.size(); pair_index++)
184 template<
unsigned DIM>
190 c_vector<double, DIM+1> phi;
192 double total_electrode_flux = 0.0;
201 unsigned node_index = node_iter->GetIndex();
203 double contribution_of_this_node = 0.0;
205 for (std::set<unsigned>::iterator iter = this->
mpMesh->
GetNode(node_index)->rGetContainingElementIndices().begin();
206 iter != this->
mpMesh->
GetNode(node_index)->rGetContainingElementIndices().end();
212 c_matrix<double, DIM, DIM> jacobian;
213 c_matrix<double, DIM, DIM> inverse_jacobian;
214 double jacobian_determinant;
217 double contribution_of_this_element = 0.0;
219 for (
unsigned quad_index=0; quad_index < pQuadRule->
GetNumQuadPoints(); quad_index++)
224 double interpolated_stimulus = 0.0;
226 for (
unsigned node_index_in_element = 0; node_index_in_element < p_element->
GetNumNodes(); node_index_in_element++)
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);
235 contribution_of_this_node += contribution_of_this_element;
238 total_electrode_flux += contribution_of_this_node;
242 int mpi_ret = MPI_Allreduce(&total_electrode_flux, &ret, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
243 assert(mpi_ret == MPI_SUCCESS);
245 MPI_Allreduce(&total_electrode_flux, &ret, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
251 assert(ret < DBL_MAX);
std::vector< double > & mrDurations
virtual DistributedVectorFactory * GetDistributedVectorFactory()
std::vector< double > mMagnitudesElectrode2
void SetCompatibleExtracellularStimulus()
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)
Node< SPACE_DIM > * GetNode(unsigned index) const
#define EXCEPTION(message)
std::vector< double > & mrMagnitudes
Element< ELEMENT_DIM, SPACE_DIM > * GetElement(unsigned index) const
double ComputeElectrodeTotalFlux(AbstractChasteRegion< DIM > *pRegion, double stimulusMagnitude)
NodeIterator GetNodeIteratorEnd()
std::vector< double > & mrPeriods
virtual unsigned GetNumNodes() const
std::vector< double > & mrStarts
virtual bool DoesContain(const ChastePoint< SPACE_DIM > &rPointToCheck) const =0
bool mGroundSecondElectrode
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
void GroundSecondElectrode(bool grounded)
void CheckForElectrodesIntersection()
std::vector< double > & mrEnds
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 ComputeBasisFunctions(const ChastePoint< ELEMENT_DIM > &rPoint, c_vector< double, ELEMENT_DIM+1 > &rReturnValue)
NodeIterator GetNodeIteratorBegin(bool skipDeletedNodes=true)
AbstractTetrahedralMesh< ELEMENT_DIM, ELEMENT_DIM > * mpMesh
~ElectrodesStimulusFactory()
unsigned GetIndex() const
ChastePoint< SPACE_DIM > GetPoint() const
double GetWeight(unsigned index) const
boost::shared_ptr< AbstractStimulusFunction > CreateStimulusForNode(Node< DIM > *pNode)
const ChastePoint< ELEMENT_DIM > & rGetQuadPoint(unsigned index) const