ElectrodesStimulusFactory.cpp

00001 /*
00002 
00003 Copyright (c) 2005-2015, University of Oxford.
00004 All rights reserved.
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Redistribution and use in source and binary forms, with or without
00013 modification, are permitted provided that the following conditions are met:
00014  * Redistributions of source code must retain the above copyright notice,
00015    this list of conditions and the following disclaimer.
00016  * Redistributions in binary form must reproduce the above copyright notice,
00017    this list of conditions and the following disclaimer in the documentation
00018    and/or other materials provided with the distribution.
00019  * Neither the name of the University of Oxford nor the names of its
00020    contributors may be used to endorse or promote products derived from this
00021    software without specific prior written permission.
00022 
00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
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); //Electrodes intersect
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)//no grounding of second electrode
00141         {
00142             //compute the two fluxes for this pair
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             //Scale the magnitude of the second electrode for this pair
00147             mMagnitudesElectrode2[pair_index] = -mMagnitudesElectrode1[pair_index]*flux_electrode_1/flux_electrode_2;
00148 
00149             // Paranoia.
00150             assert( flux_electrode_2 != 0.0);
00151             assert( ! std::isnan(mMagnitudesElectrode2[pair_index]));
00152 
00153         }
00154         else//second electrode is grounded
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//no stimulus here
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     //the basis functions
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             //loop over the elements where this node is contained
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                 /*Determine jacobian for this element*/
00213                 c_matrix<double, DIM, DIM> jacobian;
00214                 c_matrix<double, DIM, DIM> inverse_jacobian;//unused here, but needed for the function below
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;//...to this node
00219                  // loop over Gauss points
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                     //loop over nodes in this element
00227                     for (unsigned node_index_in_element = 0; node_index_in_element < p_element->GetNumNodes(); node_index_in_element++)
00228                     {
00229                         //const Node<DIM>* p_node = p_element->GetNode(node_index_in_element);
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                 }/*end of loop over gauss points*/
00236                 contribution_of_this_node += contribution_of_this_element;
00237 
00238             }/*end of loop over elements where the node is contained*/
00239             total_electrode_flux += contribution_of_this_node;
00240         }/* end of if that checks if node is in the electrode*/
00241     }/* end of loop over nodes in the mesh*/
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     //clear up memory
00250     delete pQuadRule;
00251 
00252     assert(ret < DBL_MAX);
00253     return ret;
00254 }
00255 
00257 // Explicit instantiation
00259 
00260 template class ElectrodesStimulusFactory<1>;
00261 template class ElectrodesStimulusFactory<2>;
00262 template class ElectrodesStimulusFactory<3>;
00263 

Generated by  doxygen 1.6.2