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