ElectrodesStimulusFactory.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2011
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 
00029 #include "ElectrodesStimulusFactory.hpp"
00030 #include "DistributedTetrahedralMesh.hpp"
00031 #include "IsNan.hpp"
00032 #include "HeartConfig.hpp"
00033 #include "GaussianQuadratureRule.hpp"
00034 
00035 template<unsigned DIM>
00036 ElectrodesStimulusFactory<DIM>::ElectrodesStimulusFactory(std::vector<std::pair<AbstractChasteRegion<DIM>*, AbstractChasteRegion<DIM>*> >& rElectrodePairs,
00037                                                           std::vector<double>& rStimulusMagnitudes,
00038                                                           std::vector<double>& rDurations,
00039                                                           std::vector<double>& rPeriods,
00040                                                           std::vector<double>& rStarts,
00041                                                           std::vector<double>& rEnds)
00042     : mrElectrodePairs(rElectrodePairs),
00043       mrMagnitudes(rStimulusMagnitudes),
00044       mrDurations(rDurations),
00045       mrPeriods(rPeriods),
00046       mrStarts(rStarts),
00047       mrEnds(rEnds),
00048       mGroundSecondElectrode(false)
00049 {
00050 
00051     if ( ( rElectrodePairs.size() != rStimulusMagnitudes.size() ) ||
00052        ( rElectrodePairs.size() != rDurations.size() ) ||
00053        ( rElectrodePairs.size() != rPeriods.size() ) ||
00054        ( rElectrodePairs.size() != rStarts.size() ) ||
00055        ( rElectrodePairs.size() != rEnds.size() ) )
00056     {
00057         EXCEPTION ("Vector of electrode pairs and vector of stimulation paremeters must have the same size");
00058     }
00059 
00060     mMagnitudesElectrode1 = mrMagnitudes;
00061     mMagnitudesElectrode2 = mrMagnitudes;
00062 }
00063 
00064 template<unsigned DIM>
00065 ElectrodesStimulusFactory<DIM>::~ElectrodesStimulusFactory()
00066 {
00067 }
00068 
00069 template<unsigned DIM>
00070 void ElectrodesStimulusFactory<DIM>::CheckForElectrodesIntersection()
00071 {
00072     std::vector<unsigned> nodes_in_all_electrodes;
00073     for (unsigned global_node_index = 0; global_node_index < this->mpMesh->GetNumNodes(); global_node_index++)
00074     {
00075         if (this->mpMesh->GetDistributedVectorFactory()->IsGlobalIndexLocal(global_node_index))
00076         {
00077             for (unsigned pair_index = 0; pair_index <mrElectrodePairs.size(); pair_index++)
00078             {
00079                 if ( mrElectrodePairs[pair_index].first->DoesContain( this->mpMesh->GetNode(global_node_index)->GetPoint() ) )
00080                 {
00081                     nodes_in_all_electrodes.push_back( this->mpMesh->GetNode(global_node_index)->GetIndex() );
00082                 }
00083                 if ( mrElectrodePairs[pair_index].second->DoesContain( this->mpMesh->GetNode(global_node_index)->GetPoint() ) )
00084                 {
00085                     nodes_in_all_electrodes.push_back( this->mpMesh->GetNode(global_node_index)->GetIndex() );
00086                 }
00087             }
00088         }
00089     }
00090     PetscTools::Barrier();
00091     for (unsigned node_index = 0; node_index < nodes_in_all_electrodes.size(); node_index++)
00092     {
00093         unsigned number_of_hits = 0;
00094         for (unsigned node_to_check = 0; node_to_check < nodes_in_all_electrodes.size(); node_to_check++)
00095         {
00096             if (nodes_in_all_electrodes[node_index] == nodes_in_all_electrodes[node_to_check] )
00097             {
00098                 number_of_hits++;
00099             }
00100         }
00101         if (number_of_hits>1)
00102         {
00103             EXCEPTION("Two or more electrodes intersect with each other");
00104         }
00105     }
00106 }
00107 
00108 template<unsigned DIM>
00109 void ElectrodesStimulusFactory<DIM>::GroundSecondElectrode(bool grounded)
00110 {
00111     mGroundSecondElectrode = grounded;
00112 }
00113 
00114 template<unsigned DIM>
00115 void ElectrodesStimulusFactory<DIM>::SetCompatibleExtracellularStimulus()
00116 {
00117     assert(this->mpMesh!=NULL);
00118     try
00119     {
00120         CheckForElectrodesIntersection();
00121     }
00122     catch (Exception &e)
00123     {
00124         PetscTools::ReplicateException(true); //Electrodes intersect
00125         throw e;
00126     }
00127     PetscTools::ReplicateException(false);
00128 
00129     for (unsigned pair_index = 0; pair_index < mrElectrodePairs.size(); pair_index++)
00130     {
00131 
00132         if (!mGroundSecondElectrode)//no grounding of second electrode
00133         {
00134             //compute the two fluxes for this pair
00135             double flux_electrode_1  = ComputeElectrodeTotalFlux(mrElectrodePairs[pair_index].first, mMagnitudesElectrode1[pair_index]);
00136             double flux_electrode_2  = ComputeElectrodeTotalFlux(mrElectrodePairs[pair_index].second, mMagnitudesElectrode2[pair_index]);
00137 
00138             //Scale the magnitude of the second electrode for this pair
00139             mMagnitudesElectrode2[pair_index] = -mMagnitudesElectrode1[pair_index]*flux_electrode_1/flux_electrode_2;
00140 
00141             // Paranoia.
00142             assert( flux_electrode_2 != 0.0);
00143             assert( ! std::isnan(mMagnitudesElectrode2[pair_index]));
00144 
00145         }
00146         else//second electrode is grounded
00147         {
00148             this->mGroundedRegions.push_back( mrElectrodePairs[pair_index].second );
00149         }
00150     }
00151 }
00152 
00153 template<unsigned DIM>
00154 boost::shared_ptr<AbstractStimulusFunction> ElectrodesStimulusFactory<DIM>::CreateStimulusForNode(unsigned nodeIndex)
00155 {
00156     boost::shared_ptr<RegularStimulus> p_stimulus;
00157     for (unsigned pair_index = 0; pair_index < mrElectrodePairs.size(); pair_index++)
00158     {
00159         if (mrElectrodePairs[pair_index].first->DoesContain(this->mpMesh->GetNode(nodeIndex)->GetPoint()) )
00160         {
00161             p_stimulus.reset( new RegularStimulus(mMagnitudesElectrode1[pair_index], mrDurations[pair_index], mrPeriods[pair_index], mrStarts[pair_index], mrEnds[pair_index]));
00162 
00163         }
00164         else if (mrElectrodePairs[pair_index].second->DoesContain(this->mpMesh->GetNode(nodeIndex)->GetPoint()) )
00165         {
00166             p_stimulus.reset ( new RegularStimulus(mMagnitudesElectrode2[pair_index], mrDurations[pair_index], mrPeriods[pair_index], mrStarts[pair_index], mrEnds[pair_index]));
00167 
00168         }
00169         else//no stimulus here
00170         {
00171             p_stimulus.reset ( new RegularStimulus(0.0, mrDurations[pair_index], mrPeriods[pair_index], mrStarts[pair_index], mrEnds[pair_index]));
00172         }
00173     }
00174     return p_stimulus;
00175 }
00176 
00177 template<unsigned DIM>
00178 double ElectrodesStimulusFactory<DIM>::ComputeElectrodeTotalFlux(AbstractChasteRegion<DIM>* pRegion, double stimulusMagnitude)
00179 {
00180         //assume 2 gauss points
00181         GaussianQuadratureRule<DIM>* pQuadRule = new GaussianQuadratureRule<DIM>(2u);
00182 
00183         //the basis functions
00184         c_vector<double, DIM+1> phi;
00185 
00186         double total_electrode_flux = 0.0;
00187         double ret;
00188 
00189         for (typename AbstractTetrahedralMesh<DIM,DIM>::NodeIterator iter=this->mpMesh->GetNodeIteratorBegin();
00190              iter != this->mpMesh->GetNodeIteratorEnd();
00191              ++iter)
00192         {
00193             if ( pRegion->DoesContain( (*iter).GetPoint() ) )
00194             {
00195                 unsigned node_index =     iter->GetIndex();
00196                 assert(node_index < this->mpMesh->GetNumNodes());
00197                 double contribution_of_this_node = 0.0;
00198                 //loop over the elements where this node is contained
00199                 for(std::set<unsigned>::iterator iter = this->mpMesh->GetNode(node_index)->rGetContainingElementIndices().begin();
00200                     iter != this->mpMesh->GetNode(node_index)->rGetContainingElementIndices().end();
00201                     ++iter)
00202                 {
00203                     Element<DIM,DIM>* p_element = this->mpMesh->GetElement(*iter);
00204 
00205                     /*Determine jacobian for this element*/
00206                     c_matrix<double, DIM, DIM> jacobian;
00207                     c_matrix<double, DIM, DIM> inverse_jacobian;//unused here, but needed for the function below
00208                     double jacobian_determinant;
00209                     this->mpMesh->GetInverseJacobianForElement(p_element->GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00210 
00211                     double contribution_of_this_element = 0.0;//...to this node
00212                      // loop over Gauss points
00213                     for (unsigned quad_index=0; quad_index < pQuadRule->GetNumQuadPoints(); quad_index++)
00214                     {
00215                         const ChastePoint<DIM>& quad_point = pQuadRule->rGetQuadPoint(quad_index);
00216                         BasisFunction::ComputeBasisFunctions(quad_point, phi);
00217 
00218                         double interpolated_stimulus = 0.0;
00219                         //loop over nodes in this element
00220                         for (unsigned node_index_in_element = 0; node_index_in_element < p_element->GetNumNodes(); node_index_in_element++)
00221                         {
00222                             //const Node<DIM>* p_node = p_element->GetNode(node_index_in_element);
00223                             assert(p_element->GetNumNodes() == DIM+1);
00224                             interpolated_stimulus += stimulusMagnitude*phi(node_index_in_element);
00225                             contribution_of_this_element += interpolated_stimulus * phi(node_index_in_element) * jacobian_determinant * pQuadRule->GetWeight(quad_index);
00226                         }
00227 
00228                     }/*end of loop over gauss points*/
00229                     contribution_of_this_node += contribution_of_this_element;
00230 
00231                 }/*end of loop over elements where the node is contained*/
00232                 total_electrode_flux += contribution_of_this_node;
00233             }/* end of if that checks if node is in the electrode*/
00234         }/* end of loop over nodes in the mesh*/
00235 #ifndef NDEBUG
00236         int mpi_ret = MPI_Allreduce(&total_electrode_flux, &ret, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
00237         assert(mpi_ret == MPI_SUCCESS);
00238 #else
00239         MPI_Allreduce(&total_electrode_flux, &ret, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
00240 #endif
00241 
00242         //clear up memory
00243         delete pQuadRule;
00244 
00245         assert(ret < DBL_MAX);
00246         return ret;
00247 }
00248 
00250 // Explicit instantiation
00252 
00253 template class ElectrodesStimulusFactory<1>;
00254 template class ElectrodesStimulusFactory<2>;
00255 template class ElectrodesStimulusFactory<3>;
00256 
Generated on Thu Dec 22 13:00:06 2011 for Chaste by  doxygen 1.6.3