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