StreeterFibreGenerator.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 "UblasCustomFunctions.hpp"
00037 
00038 #include "StreeterFibreGenerator.hpp"
00039 
00040 #include <cmath>
00041 #include <fstream>
00042 #include <sstream>
00043 #include "OutputFileHandler.hpp"
00044 #include "Exception.hpp"
00045 //#include "HeartRegionCodes.hpp"
00046 
00047 // Add the citation for original Streeter paper.
00048 #include "Citations.hpp"
00049 static PetscBool StreeterCite = PETSC_FALSE;
00050 const char StreeterCitation[] = "@article{streeter1969fiber,\n"
00051 "  title={Fiber orientation in the canine left ventricle during diastole and systole},\n"
00052 "  author={Streeter, Daniel D and Spotnitz, Henry M and Patel, Dali P and Ross, John and Sonnenblick, Edmund H},\n"
00053 "  journal={Circulation research},\n"
00054 "  volume={24},\n"
00055 "  number={3},\n"
00056 "  pages={339--347},\n"
00057 "  year={1969},\n"
00058 "  publisher={Am Heart Assoc}\n"
00059 "}\n";
00060 
00061 template<unsigned SPACE_DIM>
00062 double StreeterFibreGenerator<SPACE_DIM>::GetAveragedThicknessLocalNode(
00063         const unsigned nodeIndex, const std::vector<double>& wallThickness) const
00064 {
00065     if (nodeIndex < this->mpMesh->GetDistributedVectorFactory()->GetLow() ||
00066         nodeIndex >= this->mpMesh->GetDistributedVectorFactory()->GetHigh() )
00067     {
00068         return 0.0;  //Don't calculate this for nodes which aren't local
00069     }
00070 
00071     // Initialise the average with the value corresponding to the current node
00072     double average = wallThickness[nodeIndex];
00073     unsigned nodes_visited = 1;
00074 
00075     // Use a set to store visited nodes
00076     std::set<unsigned> visited_nodes;
00077     visited_nodes.insert(nodeIndex);
00078 
00079     Node<SPACE_DIM>* p_current_node = this->mpMesh->GetNode(nodeIndex);
00080 
00081     // Loop over the elements containing the given node
00082     for(typename Node<SPACE_DIM>::ContainingElementIterator element_iterator = p_current_node->ContainingElementsBegin();
00083         element_iterator != p_current_node->ContainingElementsEnd();
00084         ++element_iterator)
00085     {
00086         // Get a pointer to the container element
00087         Element<SPACE_DIM,SPACE_DIM>* p_containing_element = this->mpMesh->GetElement(*element_iterator);
00088 
00089        // Loop over the nodes of the element
00090        for(unsigned node_local_index=0;
00091            node_local_index<p_containing_element->GetNumNodes();
00092            node_local_index++)
00093        {
00094             Node<SPACE_DIM>* p_neighbour_node = p_containing_element->GetNode(node_local_index);
00095             unsigned neighbour_node_index = p_neighbour_node->GetIndex();
00096 
00097             // Check if the neighbour node has already been visited
00098             if (visited_nodes.find(neighbour_node_index) == visited_nodes.end())
00099             {
00100                 average += wallThickness[neighbour_node_index];
00101                 visited_nodes.insert(neighbour_node_index);
00102                 nodes_visited++;
00103             }
00104        }
00105     }
00106 
00107     return average/nodes_visited;
00108 }
00109 
00110 
00111 
00112 
00113 template<unsigned SPACE_DIM>
00114 double StreeterFibreGenerator<SPACE_DIM>::GetFibreMaxAngle(
00115         const c_vector<HeartRegionType, SPACE_DIM+1>& nodesRegionsForElement) const
00116 {
00117     unsigned lv=0, rv=0;
00118 
00119     for (unsigned index=0; index<SPACE_DIM+1; index++)
00120     {
00121         switch (nodesRegionsForElement[index])
00122         {
00123             case HeartGeometryInformation<SPACE_DIM>::LEFT_VENTRICLE_SURFACE:
00124             case HeartGeometryInformation<SPACE_DIM>::LEFT_VENTRICLE_WALL:
00125             case HeartGeometryInformation<SPACE_DIM>::LEFT_SEPTUM:
00126                 lv++;
00127                 break;
00128 
00129             case HeartGeometryInformation<SPACE_DIM>::RIGHT_VENTRICLE_SURFACE:
00130             case HeartGeometryInformation<SPACE_DIM>::RIGHT_VENTRICLE_WALL:
00131             case HeartGeometryInformation<SPACE_DIM>::RIGHT_SEPTUM:
00132                 rv++;
00133                 break;
00134 
00135             case HeartGeometryInformation<SPACE_DIM>::UNKNOWN:
00136             default:
00137                 NEVER_REACHED;
00138         }
00139     }
00140 
00141     // If most of the nodes are in the right ventricle
00142     if (rv>lv)
00143     {
00144         return M_PI/4.0;
00145     }
00146 
00147     // Anywhere else
00148     return M_PI/3.0;
00149 }
00150 
00151 template<unsigned SPACE_DIM>
00152 StreeterFibreGenerator<SPACE_DIM>::StreeterFibreGenerator(AbstractTetrahedralMesh<SPACE_DIM,SPACE_DIM>& rMesh)
00153     : AbstractPerElementWriter<SPACE_DIM,SPACE_DIM,SPACE_DIM*SPACE_DIM>(&rMesh),
00154       mpGeometryInfo(NULL),
00155       mApexToBase(zero_vector<double>(SPACE_DIM)),
00156       mLogInfo(false)
00157 {
00158     // Record a reference for the calculations performed here, can be extracted with the '-citations' flag.
00159     Citations::Register(StreeterCitation, &StreeterCite);
00160 
00161     mWallThickness.resize(rMesh.GetNumNodes());
00162     mAveragedWallThickness.resize(rMesh.GetNumNodes());
00163 }
00164 
00165 template<unsigned SPACE_DIM>
00166 StreeterFibreGenerator<SPACE_DIM>::~StreeterFibreGenerator()
00167 {
00168     delete mpGeometryInfo;
00169 }
00170 
00171 template<unsigned SPACE_DIM>
00172 void StreeterFibreGenerator<SPACE_DIM>::SetSurfaceFiles(
00173             const std::string &rEpicardiumFile,
00174             const std::string &rRightVentricleFile,
00175             const std::string &rLeftVentricleFile,
00176             bool indexFromZero)
00177 {
00178     // Compute the distance map of each surface
00179      mpGeometryInfo = new HeartGeometryInformation<SPACE_DIM>(*(this->mpMesh), rEpicardiumFile, rLeftVentricleFile, rRightVentricleFile, indexFromZero);
00180 }
00181 
00182 template<unsigned SPACE_DIM>
00183 void StreeterFibreGenerator<SPACE_DIM>::WriteHeaderOnMaster()
00184 {
00185     *(this->mpMasterFile) << this->mpMesh->GetNumElements();
00186     *(this->mpMasterFile) << std::setprecision(16);
00187 }
00188 
00189 template<unsigned SPACE_DIM>
00190 void StreeterFibreGenerator<SPACE_DIM>::PreWriteCalculations(OutputFileHandler& rOutputDirectory)
00191 {
00192     assert(SPACE_DIM == 3);
00193     if (mpGeometryInfo == NULL)
00194     {
00195         EXCEPTION("Files defining the heart surfaces not set");
00196     }
00197 
00198     // Open files
00199     out_stream p_regions_file, p_thickness_file, p_ave_thickness_file;
00200 
00201     //Make sure that only the master process writes the log files if requested.
00202     bool logInfo = PetscTools::AmMaster() & mLogInfo;
00203 
00204     if (logInfo)
00205     {
00206         p_regions_file  = rOutputDirectory.OpenOutputFile("node_regions.data");
00207         p_thickness_file = rOutputDirectory.OpenOutputFile("wall_thickness.data");
00208         p_ave_thickness_file = rOutputDirectory.OpenOutputFile("averaged_thickness.data");
00209     }
00210 
00211     //We expect that the apex to base has been set
00212     if (fabs(norm_2(mApexToBase)) < DBL_EPSILON)
00213     {
00214         EXCEPTION("Apex to base vector has not been set");
00215     }
00216 
00217     // Compute wall thickness parameter
00218     unsigned num_nodes = this->mpMesh->GetNumNodes();
00219     for (unsigned node_index=0; node_index<num_nodes; node_index++)
00220     {
00221         double dist_epi, dist_endo;
00222 
00223         HeartRegionType node_region = mpGeometryInfo->GetHeartRegion(node_index);
00224 
00225         switch(node_region)
00226         {
00227             case HeartGeometryInformation<SPACE_DIM>::LEFT_VENTRICLE_SURFACE:
00228             case HeartGeometryInformation<SPACE_DIM>::LEFT_VENTRICLE_WALL:
00229                 dist_epi = mpGeometryInfo->rGetDistanceMapEpicardium()[node_index];
00230                 dist_endo = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
00231                 break;
00232 
00233             case HeartGeometryInformation<SPACE_DIM>::RIGHT_VENTRICLE_SURFACE:
00234             case HeartGeometryInformation<SPACE_DIM>::RIGHT_VENTRICLE_WALL:
00235                 dist_epi = mpGeometryInfo->rGetDistanceMapEpicardium()[node_index];
00236                 dist_endo = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
00237                 break;
00238 
00239             case HeartGeometryInformation<SPACE_DIM>::LEFT_SEPTUM:
00240                 dist_epi = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
00241                 dist_endo = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
00242                 break;
00243 
00244             case HeartGeometryInformation<SPACE_DIM>::RIGHT_SEPTUM:
00245                 dist_epi = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
00246                 dist_endo = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
00247                 break;
00248 
00249             case HeartGeometryInformation<SPACE_DIM>::UNKNOWN:
00250                 #define COVERAGE_IGNORE
00251                 std::cerr << "Wrong distances node: " << node_index << "\t"
00252                           << "Epi " << mpGeometryInfo->rGetDistanceMapEpicardium()[node_index] << "\t"
00253                           << "RV " << mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index] << "\t"
00254                           << "LV " << mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index]
00255                           << std::endl;
00256 
00257                 // Make wall_thickness=0 as in Martin's code
00258                 dist_epi = 1;
00259                 dist_endo = 0;
00260                 break;
00261                 #undef COVERAGE_IGNORE
00262 
00263             default:
00264                 NEVER_REACHED;
00265         }
00266 
00267         mWallThickness[node_index] = dist_endo / (dist_endo + dist_epi);
00268 
00269         if (std::isnan(mWallThickness[node_index]))
00270         {
00271             #define COVERAGE_IGNORE
00272             /*
00273              *  A node contained on both epicardium and lv (or rv) surfaces has wall thickness 0/0.
00274              *  By setting its value to 0 we consider it contained only on the lv (or rv) surface.
00275              */
00276             mWallThickness[node_index] = 0;
00277             #undef COVERAGE_IGNORE
00278         }
00279 
00280         if (logInfo)
00281         {
00282             *p_regions_file << node_region*100 << "\n";
00283             *p_thickness_file << mWallThickness[node_index] << "\n";
00284         }
00285     }
00286 
00287     /*
00288      *  For each node, average its value of e with the values of all the neighbours
00289      */
00290     std::vector<double> my_averaged_wall_thickness(num_nodes);
00291     for (unsigned node_index=0; node_index<num_nodes; node_index++)
00292     {
00293         my_averaged_wall_thickness[node_index] = GetAveragedThicknessLocalNode(node_index, mWallThickness);
00294     }
00295 
00296     // Non-local information appear as zeros in the vector
00297     MPI_Allreduce(&my_averaged_wall_thickness[0], &mAveragedWallThickness[0], num_nodes,
00298                   MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
00299 
00300     if (logInfo)
00301     {
00302         for (unsigned node_index=0; node_index<num_nodes; node_index++)
00303         {
00304              *p_ave_thickness_file << mAveragedWallThickness[node_index] << "\n";
00305         }
00306     }
00307 
00308     if (logInfo)
00309     {
00310         p_regions_file->close();
00311         p_thickness_file->close();
00312         p_ave_thickness_file->close();
00313     }
00314 }
00315 
00316 template<unsigned SPACE_DIM>
00317 void StreeterFibreGenerator<SPACE_DIM>::Visit(Element<SPACE_DIM, SPACE_DIM>* pElement,
00318                                               unsigned localElementIndex,
00319                                               c_vector<double, SPACE_DIM*SPACE_DIM>& rData)
00320 {
00321     /*
00322      *  The gradient of the averaged thickness at the element is:
00323      *
00324      *     grad_ave_wall_thickness[element_index] = ave' * BF * inv(J)
00325      *
00326      *  being : ave, averaged thickness values of the nodes defining the element
00327      *          J,   the Jacobian of the element as defined in class Element.
00328      *                               (-1 -1 -1)
00329      *          BF,  basis functions ( 1  0  0)
00330      *                               ( 0  1  0)
00331      *                               ( 0  0  1)
00332      *
00333      *  Defined as u in Streeter paper.
00334      */
00335     c_vector<double, SPACE_DIM> grad_ave_wall_thickness;
00336     c_vector<double, SPACE_DIM+1> elem_nodes_ave_thickness;
00337     double element_averaged_thickness = 0.0;
00338     c_vector<HeartRegionType, SPACE_DIM+1> elem_nodes_region;
00339 
00340     for (unsigned local_node_index=0; local_node_index<SPACE_DIM+1; local_node_index++)
00341     {
00342         // Get node's global index
00343         unsigned global_node_index = pElement->GetNode(local_node_index)->GetIndex();
00344 
00345         elem_nodes_ave_thickness[local_node_index] = mAveragedWallThickness[global_node_index];
00346         elem_nodes_region[local_node_index] = mpGeometryInfo->GetHeartRegion(global_node_index);
00347 
00348         // Calculate wall thickness averaged value for the element
00349         element_averaged_thickness +=  mWallThickness[global_node_index];
00350     }
00351 
00352     element_averaged_thickness /= SPACE_DIM+1;
00353 
00354     c_matrix<double, SPACE_DIM+1, SPACE_DIM> basis_functions( zero_matrix<double>(4u,3u) );
00355     basis_functions(0,0) = basis_functions(0,1) = basis_functions(0,2) = -1.0;
00356     basis_functions(1,0) = basis_functions(2,1) = basis_functions(3,2) =  1.0;
00357 
00358     c_matrix<double, SPACE_DIM+1, SPACE_DIM> temp;
00359     c_matrix<double, SPACE_DIM, SPACE_DIM> jacobian, inverse_jacobian;
00360     double jacobian_det;
00361     unsigned element_index = pElement->GetIndex();
00362     this->mpMesh->GetInverseJacobianForElement(element_index, jacobian, jacobian_det, inverse_jacobian);
00363     noalias(temp) = prod (basis_functions, inverse_jacobian);
00364     noalias(grad_ave_wall_thickness) = prod(elem_nodes_ave_thickness, temp);
00365 
00366     grad_ave_wall_thickness /= norm_2(grad_ave_wall_thickness);
00367 
00368     /*
00369      * Normal to the gradient (v in Streeter paper) which is then the circumferential direction
00370      * (it will be the fibre direction after rotation)
00371      *
00372      * Computed as the cross product with the base-apex direction (originally assumed base-apex axis is x). The output vector is not normal,
00373      * since the angle between them may be != 90, normalise it.
00374      */
00375     c_vector<double, SPACE_DIM> fibre_direction = VectorProduct(grad_ave_wall_thickness, mApexToBase);
00376     fibre_direction /= norm_2(fibre_direction);
00377 
00378     /*
00379      *  Longitude direction (w in Streeter paper)
00380      */
00381     c_vector<double, SPACE_DIM> longitude_direction = VectorProduct(grad_ave_wall_thickness, fibre_direction);
00382 
00383     /*
00384      *  Compute fibre to v angle: alpha = R*(1-2e)^3
00385      *
00386      *    R is the maximum angle between the fibre and the v axis (heart region dependant)
00387      *    (1 - 2e)^3 scales it by a value in [-1, 1] defining the rotation of the fibre based
00388      *       on the position in the wall
00389      */
00390     double alpha = GetFibreMaxAngle(elem_nodes_region) * SmallPow( (1 - 2*element_averaged_thickness), 3 );
00391 
00392     /*
00393      *  Apply alpha rotation about the u axis to the orthonormal basis
00394      *
00395      *               ( u(1) v(1) w(1) )
00396      *   (u, v, w) = ( u(2) v(2) w(2) )
00397      *               ( u(3) v(3) w(3) )
00398      *
00399      *  The following matrix defines a rotation about the u axis
00400      *
00401      *                 ( 1        0           0      ) (u')
00402      *   R = (u, v, w) ( 0    cos(alpha) -sin(alpha) ) (v')
00403      *                 ( 0    sin(alpha)  cos(alpha) ) (w')
00404      *
00405      *  The rotated basis is computed like:
00406      *
00407      *                                             ( 1        0           0      )
00408      *  (u, v_r, w_r ) = R * (u, v, w) = (u, v, w) ( 0    cos(alpha) -sin(alpha) )
00409      *                                             ( 0    sin(alpha)  cos(alpha) )
00410      *
00411      *  Which simplifies to:
00412      *
00413      *   v_r =  v*cos(alpha) + w*sin(alpha)
00414      *   w_r = -v*sin(alpha) + w*cos(alpha)
00415      */
00416     c_vector<double, SPACE_DIM> rotated_fibre_direction = fibre_direction*cos(alpha) + longitude_direction*sin(alpha);
00417     c_vector<double, SPACE_DIM> rotated_longitude_direction = -fibre_direction*sin(alpha) + longitude_direction*cos(alpha);
00418 
00419 
00420     /*
00421      * Test the orthonormality of the basis
00422      */
00423     assert( fabs(norm_2(rotated_fibre_direction) - 1) < 100*DBL_EPSILON );
00424     assert( fabs(norm_2(grad_ave_wall_thickness) - 1) < 100*DBL_EPSILON );
00425     assert( fabs(norm_2(rotated_longitude_direction) - 1) < 100*DBL_EPSILON );
00426 
00427     assert( fabs(inner_prod(rotated_fibre_direction, grad_ave_wall_thickness)) < 100*DBL_EPSILON );
00428     assert( fabs(inner_prod(rotated_fibre_direction, rotated_longitude_direction)) < 100*DBL_EPSILON);
00429     assert( fabs(inner_prod(grad_ave_wall_thickness, rotated_longitude_direction)) < 100*DBL_EPSILON);
00430 
00431     /*
00432      *  Output the direction of the myofibre, the transverse to it in the plane
00433      *  of the myocite laminae and the normal to this laminae (in that order)
00434      *
00435      *  See Fig. 1 "Laminar Structure of the Heart: a mathematical model" LeGrice et al. 97
00436      *
00437      */
00438     rData[0] = rotated_fibre_direction[0];
00439     rData[1] = rotated_fibre_direction[1];
00440     rData[2] = rotated_fibre_direction[2];
00441     rData[3] = grad_ave_wall_thickness[0];
00442     rData[4] = grad_ave_wall_thickness[1];
00443     rData[5] = grad_ave_wall_thickness[2];
00444     rData[6] = rotated_longitude_direction[0];
00445     rData[7] = rotated_longitude_direction[1];
00446     rData[8] = rotated_longitude_direction[2];
00447 }
00448 
00449 
00450 template<unsigned SPACE_DIM>
00451 void StreeterFibreGenerator<SPACE_DIM>::SetApexToBase(const c_vector<double, SPACE_DIM>& apexToBase)
00452 {
00453     double norm = norm_2(apexToBase);
00454     if (norm < DBL_EPSILON)
00455     {
00456         EXCEPTION("Apex to base vector should be non-zero");
00457     }
00458     mApexToBase = apexToBase / norm;
00459 }
00460 
00461 template<unsigned SPACE_DIM>
00462 void StreeterFibreGenerator<SPACE_DIM>::SetApexToBase(unsigned axis)
00463 {
00464     if (axis >= SPACE_DIM)
00465     {
00466         EXCEPTION("Apex to base coordinate axis was out of range");
00467     }
00468     mApexToBase = zero_vector<double>(SPACE_DIM);
00469     mApexToBase[axis] = 1.0;
00470 }
00471 
00472 template<unsigned SPACE_DIM>
00473 void StreeterFibreGenerator<SPACE_DIM>::SetLogInfo(bool logInfo)
00474 {
00475     mLogInfo = logInfo;
00476 }
00477 
00479 // Explicit instantiation
00481 #define COVERAGE_IGNORE
00482 template class StreeterFibreGenerator<3>;
00483 #undef COVERAGE_IGNORE

Generated by  doxygen 1.6.2