00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036 #include "ElectrodesStimulusFactory.hpp"
00037 #include "DistributedTetrahedralMesh.hpp"
00038 #include "IsNan.hpp"
00039 #include "HeartConfig.hpp"
00040 #include "GaussianQuadratureRule.hpp"
00041
00042 template<unsigned DIM>
00043 ElectrodesStimulusFactory<DIM>::ElectrodesStimulusFactory(std::vector<std::pair<AbstractChasteRegion<DIM>*, AbstractChasteRegion<DIM>*> >& rElectrodePairs,
00044 std::vector<double>& rStimulusMagnitudes,
00045 std::vector<double>& rDurations,
00046 std::vector<double>& rPeriods,
00047 std::vector<double>& rStarts,
00048 std::vector<double>& rEnds)
00049 : mrElectrodePairs(rElectrodePairs),
00050 mrMagnitudes(rStimulusMagnitudes),
00051 mrDurations(rDurations),
00052 mrPeriods(rPeriods),
00053 mrStarts(rStarts),
00054 mrEnds(rEnds),
00055 mGroundSecondElectrode(false)
00056 {
00057
00058 if ( ( rElectrodePairs.size() != rStimulusMagnitudes.size() ) ||
00059 ( rElectrodePairs.size() != rDurations.size() ) ||
00060 ( rElectrodePairs.size() != rPeriods.size() ) ||
00061 ( rElectrodePairs.size() != rStarts.size() ) ||
00062 ( rElectrodePairs.size() != rEnds.size() ) )
00063 {
00064 EXCEPTION ("Vector of electrode pairs and vector of stimulation paremeters must have the same size");
00065 }
00066
00067 mMagnitudesElectrode1 = mrMagnitudes;
00068 mMagnitudesElectrode2 = mrMagnitudes;
00069 }
00070
00071 template<unsigned DIM>
00072 ElectrodesStimulusFactory<DIM>::~ElectrodesStimulusFactory()
00073 {
00074 }
00075
00076 template<unsigned DIM>
00077 void ElectrodesStimulusFactory<DIM>::CheckForElectrodesIntersection()
00078 {
00079 std::vector<unsigned> nodes_in_all_electrodes;
00080 for (unsigned global_node_index = 0; global_node_index < this->mpMesh->GetNumNodes(); global_node_index++)
00081 {
00082 if (this->mpMesh->GetDistributedVectorFactory()->IsGlobalIndexLocal(global_node_index))
00083 {
00084 for (unsigned pair_index = 0; pair_index <mrElectrodePairs.size(); pair_index++)
00085 {
00086 if ( mrElectrodePairs[pair_index].first->DoesContain( this->mpMesh->GetNode(global_node_index)->GetPoint() ) )
00087 {
00088 nodes_in_all_electrodes.push_back( this->mpMesh->GetNode(global_node_index)->GetIndex() );
00089 }
00090 if ( mrElectrodePairs[pair_index].second->DoesContain( this->mpMesh->GetNode(global_node_index)->GetPoint() ) )
00091 {
00092 nodes_in_all_electrodes.push_back( this->mpMesh->GetNode(global_node_index)->GetIndex() );
00093 }
00094 }
00095 }
00096 }
00097 PetscTools::Barrier();
00098 for (unsigned node_index = 0; node_index < nodes_in_all_electrodes.size(); node_index++)
00099 {
00100 unsigned number_of_hits = 0;
00101 for (unsigned node_to_check = 0; node_to_check < nodes_in_all_electrodes.size(); node_to_check++)
00102 {
00103 if (nodes_in_all_electrodes[node_index] == nodes_in_all_electrodes[node_to_check] )
00104 {
00105 number_of_hits++;
00106 }
00107 }
00108 if (number_of_hits>1)
00109 {
00110 EXCEPTION("Two or more electrodes intersect with each other");
00111 }
00112 }
00113 }
00114
00115 template<unsigned DIM>
00116 void ElectrodesStimulusFactory<DIM>::GroundSecondElectrode(bool grounded)
00117 {
00118 mGroundSecondElectrode = grounded;
00119 }
00120
00121 template<unsigned DIM>
00122 void ElectrodesStimulusFactory<DIM>::SetCompatibleExtracellularStimulus()
00123 {
00124 assert(this->mpMesh!=NULL);
00125 try
00126 {
00127 CheckForElectrodesIntersection();
00128 }
00129 catch (Exception &e)
00130 {
00131 PetscTools::ReplicateException(true);
00132 throw e;
00133 }
00134 PetscTools::ReplicateException(false);
00135
00136 for (unsigned pair_index = 0; pair_index < mrElectrodePairs.size(); pair_index++)
00137 {
00138
00139 if (!mGroundSecondElectrode)
00140 {
00141
00142 double flux_electrode_1 = ComputeElectrodeTotalFlux(mrElectrodePairs[pair_index].first, mMagnitudesElectrode1[pair_index]);
00143 double flux_electrode_2 = ComputeElectrodeTotalFlux(mrElectrodePairs[pair_index].second, mMagnitudesElectrode2[pair_index]);
00144
00145
00146 mMagnitudesElectrode2[pair_index] = -mMagnitudesElectrode1[pair_index]*flux_electrode_1/flux_electrode_2;
00147
00148
00149 assert( flux_electrode_2 != 0.0);
00150 assert( ! std::isnan(mMagnitudesElectrode2[pair_index]));
00151
00152 }
00153 else
00154 {
00155 this->mGroundedRegions.push_back( mrElectrodePairs[pair_index].second );
00156 }
00157 }
00158 }
00159
00160 template<unsigned DIM>
00161 boost::shared_ptr<AbstractStimulusFunction> ElectrodesStimulusFactory<DIM>::CreateStimulusForNode(Node<DIM>* pNode)
00162 {
00163 boost::shared_ptr<RegularStimulus> p_stimulus;
00164 for (unsigned pair_index = 0; pair_index < mrElectrodePairs.size(); pair_index++)
00165 {
00166 if (mrElectrodePairs[pair_index].first->DoesContain(pNode->GetPoint()) )
00167 {
00168 p_stimulus.reset( new RegularStimulus(mMagnitudesElectrode1[pair_index], mrDurations[pair_index], mrPeriods[pair_index], mrStarts[pair_index], mrEnds[pair_index]));
00169
00170 }
00171 else if (mrElectrodePairs[pair_index].second->DoesContain(pNode->GetPoint()) )
00172 {
00173 p_stimulus.reset ( new RegularStimulus(mMagnitudesElectrode2[pair_index], mrDurations[pair_index], mrPeriods[pair_index], mrStarts[pair_index], mrEnds[pair_index]));
00174
00175 }
00176 else
00177 {
00178 p_stimulus.reset ( new RegularStimulus(0.0, mrDurations[pair_index], mrPeriods[pair_index], mrStarts[pair_index], mrEnds[pair_index]));
00179 }
00180 }
00181 return p_stimulus;
00182 }
00183
00184 template<unsigned DIM>
00185 double ElectrodesStimulusFactory<DIM>::ComputeElectrodeTotalFlux(AbstractChasteRegion<DIM>* pRegion, double stimulusMagnitude)
00186 {
00187 GaussianQuadratureRule<DIM>* pQuadRule = new GaussianQuadratureRule<DIM>(2);
00188
00189
00190 c_vector<double, DIM+1> phi;
00191
00192 double total_electrode_flux = 0.0;
00193 double ret;
00194
00195 for (typename AbstractTetrahedralMesh<DIM,DIM>::NodeIterator node_iter=this->mpMesh->GetNodeIteratorBegin();
00196 node_iter != this->mpMesh->GetNodeIteratorEnd();
00197 ++node_iter)
00198 {
00199 if ( pRegion->DoesContain( (*node_iter).GetPoint() ) )
00200 {
00201 unsigned node_index = node_iter->GetIndex();
00202 assert(node_index < this->mpMesh->GetNumNodes());
00203 double contribution_of_this_node = 0.0;
00204
00205 for(std::set<unsigned>::iterator iter = this->mpMesh->GetNode(node_index)->rGetContainingElementIndices().begin();
00206 iter != this->mpMesh->GetNode(node_index)->rGetContainingElementIndices().end();
00207 ++iter)
00208 {
00209 Element<DIM,DIM>* p_element = this->mpMesh->GetElement(*iter);
00210
00211
00212 c_matrix<double, DIM, DIM> jacobian;
00213 c_matrix<double, DIM, DIM> inverse_jacobian;
00214 double jacobian_determinant;
00215 this->mpMesh->GetInverseJacobianForElement(p_element->GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00216
00217 double contribution_of_this_element = 0.0;
00218
00219 for (unsigned quad_index=0; quad_index < pQuadRule->GetNumQuadPoints(); quad_index++)
00220 {
00221 const ChastePoint<DIM>& quad_point = pQuadRule->rGetQuadPoint(quad_index);
00222 BasisFunction::ComputeBasisFunctions(quad_point, phi);
00223
00224 double interpolated_stimulus = 0.0;
00225
00226 for (unsigned node_index_in_element = 0; node_index_in_element < p_element->GetNumNodes(); node_index_in_element++)
00227 {
00228
00229 assert(p_element->GetNumNodes() == DIM+1);
00230 interpolated_stimulus += stimulusMagnitude*phi(node_index_in_element);
00231 contribution_of_this_element += interpolated_stimulus * phi(node_index_in_element) * jacobian_determinant * pQuadRule->GetWeight(quad_index);
00232 }
00233
00234 }
00235 contribution_of_this_node += contribution_of_this_element;
00236
00237 }
00238 total_electrode_flux += contribution_of_this_node;
00239 }
00240 }
00241 #ifndef NDEBUG
00242 int mpi_ret = MPI_Allreduce(&total_electrode_flux, &ret, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
00243 assert(mpi_ret == MPI_SUCCESS);
00244 #else
00245 MPI_Allreduce(&total_electrode_flux, &ret, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
00246 #endif
00247
00248
00249 delete pQuadRule;
00250
00251 assert(ret < DBL_MAX);
00252 return ret;
00253 }
00254
00256
00258
00259 template class ElectrodesStimulusFactory<1>;
00260 template class ElectrodesStimulusFactory<2>;
00261 template class ElectrodesStimulusFactory<3>;
00262