Chaste  Release::3.4
StreeterFibreGenerator.cpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34 */
35 
36 #include "UblasCustomFunctions.hpp"
37 
38 #include "StreeterFibreGenerator.hpp"
39 
40 #include <cmath>
41 #include <fstream>
42 #include <sstream>
43 #include "OutputFileHandler.hpp"
44 #include "Exception.hpp"
45 //#include "HeartRegionCodes.hpp"
46 
47 // Add the citation for original Streeter paper.
48 #include "Citations.hpp"
49 static PetscBool StreeterCite = PETSC_FALSE;
50 const char StreeterCitation[] = "@article{streeter1969fiber,\n"
51 " title={Fiber orientation in the canine left ventricle during diastole and systole},\n"
52 " author={Streeter, Daniel D and Spotnitz, Henry M and Patel, Dali P and Ross, John and Sonnenblick, Edmund H},\n"
53 " journal={Circulation research},\n"
54 " volume={24},\n"
55 " number={3},\n"
56 " pages={339--347},\n"
57 " year={1969},\n"
58 " publisher={Am Heart Assoc}\n"
59 "}\n";
60 
61 template<unsigned SPACE_DIM>
63  const unsigned nodeIndex, const std::vector<double>& wallThickness) const
64 {
65  if (nodeIndex < this->mpMesh->GetDistributedVectorFactory()->GetLow() ||
66  nodeIndex >= this->mpMesh->GetDistributedVectorFactory()->GetHigh() )
67  {
68  return 0.0; //Don't calculate this for nodes which aren't local
69  }
70 
71  // Initialise the average with the value corresponding to the current node
72  double average = wallThickness[nodeIndex];
73  unsigned nodes_visited = 1;
74 
75  // Use a set to store visited nodes
76  std::set<unsigned> visited_nodes;
77  visited_nodes.insert(nodeIndex);
78 
79  Node<SPACE_DIM>* p_current_node = this->mpMesh->GetNode(nodeIndex);
80 
81  // Loop over the elements containing the given node
82  for(typename Node<SPACE_DIM>::ContainingElementIterator element_iterator = p_current_node->ContainingElementsBegin();
83  element_iterator != p_current_node->ContainingElementsEnd();
84  ++element_iterator)
85  {
86  // Get a pointer to the container element
87  Element<SPACE_DIM,SPACE_DIM>* p_containing_element = this->mpMesh->GetElement(*element_iterator);
88 
89  // Loop over the nodes of the element
90  for(unsigned node_local_index=0;
91  node_local_index<p_containing_element->GetNumNodes();
92  node_local_index++)
93  {
94  Node<SPACE_DIM>* p_neighbour_node = p_containing_element->GetNode(node_local_index);
95  unsigned neighbour_node_index = p_neighbour_node->GetIndex();
96 
97  // Check if the neighbour node has already been visited
98  if (visited_nodes.find(neighbour_node_index) == visited_nodes.end())
99  {
100  average += wallThickness[neighbour_node_index];
101  visited_nodes.insert(neighbour_node_index);
102  nodes_visited++;
103  }
104  }
105  }
106 
107  return average/nodes_visited;
108 }
109 
110 
111 
112 
113 template<unsigned SPACE_DIM>
115  const c_vector<HeartRegionType, SPACE_DIM+1>& nodesRegionsForElement) const
116 {
117  unsigned lv=0, rv=0;
118 
119  for (unsigned index=0; index<SPACE_DIM+1; index++)
120  {
121  switch (nodesRegionsForElement[index])
122  {
126  lv++;
127  break;
128 
132  rv++;
133  break;
134 
136  default:
138  }
139  }
140 
141  // If most of the nodes are in the right ventricle
142  if (rv>lv)
143  {
144  return M_PI/4.0;
145  }
146 
147  // Anywhere else
148  return M_PI/3.0;
149 }
150 
151 template<unsigned SPACE_DIM>
153  : AbstractPerElementWriter<SPACE_DIM,SPACE_DIM,SPACE_DIM*SPACE_DIM>(&rMesh),
154  mpGeometryInfo(NULL),
155  mApexToBase(zero_vector<double>(SPACE_DIM)),
156  mLogInfo(false)
157 {
158  // Record a reference for the calculations performed here, can be extracted with the '-citations' flag.
159  Citations::Register(StreeterCitation, &StreeterCite);
160 
161  mWallThickness.resize(rMesh.GetNumNodes());
162  mAveragedWallThickness.resize(rMesh.GetNumNodes());
163 }
164 
165 template<unsigned SPACE_DIM>
167 {
168  delete mpGeometryInfo;
169 }
170 
171 template<unsigned SPACE_DIM>
173  const std::string& rEpicardiumFile,
174  const std::string& rRightVentricleFile,
175  const std::string& rLeftVentricleFile,
176  bool indexFromZero)
177 {
178  // Compute the distance map of each surface
179  mpGeometryInfo = new HeartGeometryInformation<SPACE_DIM>(*(this->mpMesh), rEpicardiumFile, rLeftVentricleFile, rRightVentricleFile, indexFromZero);
180 }
181 
182 template<unsigned SPACE_DIM>
184 {
185  *(this->mpMasterFile) << this->mpMesh->GetNumElements();
186  *(this->mpMasterFile) << std::setprecision(16);
187 }
188 
189 template<unsigned SPACE_DIM>
191 {
192  assert(SPACE_DIM == 3);
193  if (mpGeometryInfo == NULL)
194  {
195  EXCEPTION("Files defining the heart surfaces not set");
196  }
197 
198  // Open files
199  out_stream p_regions_file, p_thickness_file, p_ave_thickness_file;
200 
201  //Make sure that only the master process writes the log files if requested.
202  bool logInfo = PetscTools::AmMaster() & mLogInfo;
203 
204  if (logInfo)
205  {
206  p_regions_file = rOutputDirectory.OpenOutputFile("node_regions.data");
207  p_thickness_file = rOutputDirectory.OpenOutputFile("wall_thickness.data");
208  p_ave_thickness_file = rOutputDirectory.OpenOutputFile("averaged_thickness.data");
209  }
210 
211  //We expect that the apex to base has been set
212  if (fabs(norm_2(mApexToBase)) < DBL_EPSILON)
213  {
214  EXCEPTION("Apex to base vector has not been set");
215  }
216 
217  // Compute wall thickness parameter
218  unsigned num_nodes = this->mpMesh->GetNumNodes();
219  for (unsigned node_index=0; node_index<num_nodes; node_index++)
220  {
221  double dist_epi, dist_endo;
222 
223  HeartRegionType node_region = mpGeometryInfo->GetHeartRegion(node_index);
224 
225  switch(node_region)
226  {
229  dist_epi = mpGeometryInfo->rGetDistanceMapEpicardium()[node_index];
230  dist_endo = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
231  break;
232 
235  dist_epi = mpGeometryInfo->rGetDistanceMapEpicardium()[node_index];
236  dist_endo = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
237  break;
238 
240  dist_epi = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
241  dist_endo = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
242  break;
243 
245  dist_epi = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
246  dist_endo = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
247  break;
248 
250  #define COVERAGE_IGNORE
251  std::cerr << "Wrong distances node: " << node_index << "\t"
252  << "Epi " << mpGeometryInfo->rGetDistanceMapEpicardium()[node_index] << "\t"
253  << "RV " << mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index] << "\t"
254  << "LV " << mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index]
255  << std::endl;
256 
257  // Make wall_thickness=0 as in Martin's code
258  dist_epi = 1;
259  dist_endo = 0;
260  break;
261  #undef COVERAGE_IGNORE
262 
263  default:
265  }
266 
267  mWallThickness[node_index] = dist_endo / (dist_endo + dist_epi);
268 
269  if (std::isnan(mWallThickness[node_index]))
270  {
271  #define COVERAGE_IGNORE
272  /*
273  * A node contained on both epicardium and lv (or rv) surfaces has wall thickness 0/0.
274  * By setting its value to 0 we consider it contained only on the lv (or rv) surface.
275  */
276  mWallThickness[node_index] = 0;
277  #undef COVERAGE_IGNORE
278  }
279 
280  if (logInfo)
281  {
282  *p_regions_file << node_region*100 << "\n";
283  *p_thickness_file << mWallThickness[node_index] << "\n";
284  }
285  }
286 
287  /*
288  * For each node, average its value of e with the values of all the neighbours
289  */
290  std::vector<double> my_averaged_wall_thickness(num_nodes);
291  for (unsigned node_index=0; node_index<num_nodes; node_index++)
292  {
293  my_averaged_wall_thickness[node_index] = GetAveragedThicknessLocalNode(node_index, mWallThickness);
294  }
295 
296  // Non-local information appear as zeros in the vector
297  MPI_Allreduce(&my_averaged_wall_thickness[0], &mAveragedWallThickness[0], num_nodes,
298  MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
299 
300  if (logInfo)
301  {
302  for (unsigned node_index=0; node_index<num_nodes; node_index++)
303  {
304  *p_ave_thickness_file << mAveragedWallThickness[node_index] << "\n";
305  }
306  }
307 
308  if (logInfo)
309  {
310  p_regions_file->close();
311  p_thickness_file->close();
312  p_ave_thickness_file->close();
313  }
314 }
315 
316 template<unsigned SPACE_DIM>
318  unsigned localElementIndex,
319  c_vector<double, SPACE_DIM*SPACE_DIM>& rData)
320 {
321  /*
322  * The gradient of the averaged thickness at the element is:
323  *
324  * grad_ave_wall_thickness[element_index] = ave' * BF * inv(J)
325  *
326  * being : ave, averaged thickness values of the nodes defining the element
327  * J, the Jacobian of the element as defined in class Element.
328  * (-1 -1 -1)
329  * BF, basis functions ( 1 0 0)
330  * ( 0 1 0)
331  * ( 0 0 1)
332  *
333  * Defined as u in Streeter paper.
334  */
335  c_vector<double, SPACE_DIM> grad_ave_wall_thickness;
336  c_vector<double, SPACE_DIM+1> elem_nodes_ave_thickness;
337  double element_averaged_thickness = 0.0;
338  c_vector<HeartRegionType, SPACE_DIM+1> elem_nodes_region;
339 
340  for (unsigned local_node_index=0; local_node_index<SPACE_DIM+1; local_node_index++)
341  {
342  // Get node's global index
343  unsigned global_node_index = pElement->GetNode(local_node_index)->GetIndex();
344 
345  elem_nodes_ave_thickness[local_node_index] = mAveragedWallThickness[global_node_index];
346  elem_nodes_region[local_node_index] = mpGeometryInfo->GetHeartRegion(global_node_index);
347 
348  // Calculate wall thickness averaged value for the element
349  element_averaged_thickness += mWallThickness[global_node_index];
350  }
351 
352  element_averaged_thickness /= SPACE_DIM+1;
353 
354  c_matrix<double, SPACE_DIM+1, SPACE_DIM> basis_functions( zero_matrix<double>(4u,3u) );
355  basis_functions(0,0) = basis_functions(0,1) = basis_functions(0,2) = -1.0;
356  basis_functions(1,0) = basis_functions(2,1) = basis_functions(3,2) = 1.0;
357 
358  c_matrix<double, SPACE_DIM+1, SPACE_DIM> temp;
359  c_matrix<double, SPACE_DIM, SPACE_DIM> jacobian, inverse_jacobian;
360  double jacobian_det;
361  unsigned element_index = pElement->GetIndex();
362  this->mpMesh->GetInverseJacobianForElement(element_index, jacobian, jacobian_det, inverse_jacobian);
363  noalias(temp) = prod (basis_functions, inverse_jacobian);
364  noalias(grad_ave_wall_thickness) = prod(elem_nodes_ave_thickness, temp);
365 
366  grad_ave_wall_thickness /= norm_2(grad_ave_wall_thickness);
367 
368  /*
369  * Normal to the gradient (v in Streeter paper) which is then the circumferential direction
370  * (it will be the fibre direction after rotation)
371  *
372  * Computed as the cross product with the base-apex direction (originally assumed base-apex axis is x). The output vector is not normal,
373  * since the angle between them may be != 90, normalise it.
374  */
375  c_vector<double, SPACE_DIM> fibre_direction = VectorProduct(grad_ave_wall_thickness, mApexToBase);
376  fibre_direction /= norm_2(fibre_direction);
377 
378  /*
379  * Longitude direction (w in Streeter paper)
380  */
381  c_vector<double, SPACE_DIM> longitude_direction = VectorProduct(grad_ave_wall_thickness, fibre_direction);
382 
383  /*
384  * Compute fibre to v angle: alpha = R*(1-2e)^3
385  *
386  * R is the maximum angle between the fibre and the v axis (heart region dependant)
387  * (1 - 2e)^3 scales it by a value in [-1, 1] defining the rotation of the fibre based
388  * on the position in the wall
389  */
390  double alpha = GetFibreMaxAngle(elem_nodes_region) * SmallPow( (1 - 2*element_averaged_thickness), 3 );
391 
392  /*
393  * Apply alpha rotation about the u axis to the orthonormal basis
394  *
395  * ( u(1) v(1) w(1) )
396  * (u, v, w) = ( u(2) v(2) w(2) )
397  * ( u(3) v(3) w(3) )
398  *
399  * The following matrix defines a rotation about the u axis
400  *
401  * ( 1 0 0 ) (u')
402  * R = (u, v, w) ( 0 cos(alpha) -sin(alpha) ) (v')
403  * ( 0 sin(alpha) cos(alpha) ) (w')
404  *
405  * The rotated basis is computed like:
406  *
407  * ( 1 0 0 )
408  * (u, v_r, w_r ) = R * (u, v, w) = (u, v, w) ( 0 cos(alpha) -sin(alpha) )
409  * ( 0 sin(alpha) cos(alpha) )
410  *
411  * Which simplifies to:
412  *
413  * v_r = v*cos(alpha) + w*sin(alpha)
414  * w_r = -v*sin(alpha) + w*cos(alpha)
415  */
416  c_vector<double, SPACE_DIM> rotated_fibre_direction = fibre_direction*cos(alpha) + longitude_direction*sin(alpha);
417  c_vector<double, SPACE_DIM> rotated_longitude_direction = -fibre_direction*sin(alpha) + longitude_direction*cos(alpha);
418 
419 
420  /*
421  * Test the orthonormality of the basis
422  */
423  assert( fabs(norm_2(rotated_fibre_direction) - 1) < 100*DBL_EPSILON );
424  assert( fabs(norm_2(grad_ave_wall_thickness) - 1) < 100*DBL_EPSILON );
425  assert( fabs(norm_2(rotated_longitude_direction) - 1) < 100*DBL_EPSILON );
426 
427  assert( fabs(inner_prod(rotated_fibre_direction, grad_ave_wall_thickness)) < 100*DBL_EPSILON );
428  assert( fabs(inner_prod(rotated_fibre_direction, rotated_longitude_direction)) < 100*DBL_EPSILON);
429  assert( fabs(inner_prod(grad_ave_wall_thickness, rotated_longitude_direction)) < 100*DBL_EPSILON);
430 
431  /*
432  * Output the direction of the myofibre, the transverse to it in the plane
433  * of the myocite laminae and the normal to this laminae (in that order)
434  *
435  * See Fig. 1 "Laminar Structure of the Heart: a mathematical model" LeGrice et al. 97
436  *
437  */
438  rData[0] = rotated_fibre_direction[0];
439  rData[1] = rotated_fibre_direction[1];
440  rData[2] = rotated_fibre_direction[2];
441  rData[3] = grad_ave_wall_thickness[0];
442  rData[4] = grad_ave_wall_thickness[1];
443  rData[5] = grad_ave_wall_thickness[2];
444  rData[6] = rotated_longitude_direction[0];
445  rData[7] = rotated_longitude_direction[1];
446  rData[8] = rotated_longitude_direction[2];
447 }
448 
449 
450 template<unsigned SPACE_DIM>
451 void StreeterFibreGenerator<SPACE_DIM>::SetApexToBase(const c_vector<double, SPACE_DIM>& apexToBase)
452 {
453  double norm = norm_2(apexToBase);
454  if (norm < DBL_EPSILON)
455  {
456  EXCEPTION("Apex to base vector should be non-zero");
457  }
458  mApexToBase = apexToBase / norm;
459 }
460 
461 template<unsigned SPACE_DIM>
463 {
464  if (axis >= SPACE_DIM)
465  {
466  EXCEPTION("Apex to base coordinate axis was out of range");
467  }
468  mApexToBase = zero_vector<double>(SPACE_DIM);
469  mApexToBase[axis] = 1.0;
470 }
471 
472 template<unsigned SPACE_DIM>
474 {
475  mLogInfo = logInfo;
476 }
477 
479 // Explicit instantiation
481 #define COVERAGE_IGNORE
482 template class StreeterFibreGenerator<3>;
483 #undef COVERAGE_IGNORE
void SetLogInfo(bool logInfo=true)
Definition: Node.hpp:58
void SetApexToBase(const c_vector< double, SPACE_DIM > &apexToBase)
double SmallPow(double x, unsigned exponent)
#define EXCEPTION(message)
Definition: Exception.hpp:143
static bool AmMaster()
Definition: PetscTools.cpp:120
c_vector< T, 1 > VectorProduct(const c_vector< T, 1 > &rA, const c_vector< T, 1 > &rB)
virtual unsigned GetNumNodes() const
std::vector< double > mAveragedWallThickness
#define NEVER_REACHED
Definition: Exception.hpp:206
std::vector< double > & rGetDistanceMapEpicardium()
Node< SPACE_DIM > * GetNode(unsigned localIndex) const
out_stream OpenOutputFile(const std::string &rFileName, std::ios_base::openmode mode=std::ios::out|std::ios::trunc) const
StreeterFibreGenerator(AbstractTetrahedralMesh< SPACE_DIM, SPACE_DIM > &rMesh)
ContainingElementIterator ContainingElementsBegin() const
Definition: Node.hpp:444
double GetAveragedThicknessLocalNode(const unsigned nodeIndex, const std::vector< double > &wallThickness) const
unsigned GetNumNodes() const
void SetSurfaceFiles(const std::string &rEpicardiumFile, const std::string &rRightVentricleFile, const std::string &rLeftVentricleFile, bool indexFromZero)
double GetFibreMaxAngle(const c_vector< HeartRegionType, SPACE_DIM+1 > &nodesRegionsForElement) const
void PreWriteCalculations(OutputFileHandler &rOutputDirectory)
static void Register(const char *pCitation, PetscBool *pSet)
Definition: Citations.cpp:43
std::vector< double > mWallThickness
unsigned GetIndex() const
unsigned GetIndex() const
Definition: Node.cpp:159
ContainingElementIterator ContainingElementsEnd() const
Definition: Node.hpp:452
void Visit(Element< SPACE_DIM, SPACE_DIM > *pElement, unsigned localElementIndex, c_vector< double, SPACE_DIM *SPACE_DIM > &rData)