00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #include "UblasCustomFunctions.hpp"
00030
00031 #include "StreeterFibreGenerator.hpp"
00032
00033 #include <cmath>
00034 #include <fstream>
00035 #include <sstream>
00036 #include "OutputFileHandler.hpp"
00037 #include "Exception.hpp"
00038 #include "HeartRegionCodes.hpp"
00039
00040 template<unsigned SPACE_DIM>
00041 double StreeterFibreGenerator<SPACE_DIM>::GetAveragedThickness(
00042 const unsigned nodeIndex, const std::vector<double>& wallThickness) const
00043 {
00044
00045 double average = wallThickness[nodeIndex];
00046 unsigned nodes_visited = 1;
00047
00048
00049 std::set<unsigned> visited_nodes;
00050 visited_nodes.insert(nodeIndex);
00051
00052 Node<SPACE_DIM>* p_current_node = mrMesh.GetNode(nodeIndex);
00053
00055
00056
00057 for(typename Node<SPACE_DIM>::ContainingElementIterator element_iterator = p_current_node->ContainingElementsBegin();
00058 element_iterator != p_current_node->ContainingElementsEnd();
00059 ++element_iterator)
00060 {
00061
00062 Element<SPACE_DIM,SPACE_DIM>* p_containing_element = mrMesh.GetElement(*element_iterator);
00063
00064
00065 for(unsigned node_local_index=0;
00066 node_local_index<p_containing_element->GetNumNodes();
00067 node_local_index++)
00068 {
00069 Node<SPACE_DIM>* p_neighbour_node = p_containing_element->GetNode(node_local_index);
00070 unsigned neighbour_node_index = p_neighbour_node->GetIndex();
00071
00072
00073 if (visited_nodes.find(neighbour_node_index) == visited_nodes.end())
00074 {
00075 average += wallThickness[neighbour_node_index];
00076 visited_nodes.insert(neighbour_node_index);
00077 nodes_visited++;
00078 }
00079 }
00080 }
00081
00082 return average/nodes_visited;
00083 }
00084
00085
00086
00087
00088 template<unsigned SPACE_DIM>
00089 double StreeterFibreGenerator<SPACE_DIM>::GetFibreMaxAngle(
00090 const c_vector<HeartRegionType, SPACE_DIM+1>& nodesRegionsForElement) const
00091 {
00092 unsigned lv=0, rv=0;
00093
00094 for (unsigned index=0; index<SPACE_DIM+1; index++)
00095 {
00096 switch (nodesRegionsForElement[index])
00097 {
00098 case HeartRegionCode::LEFT_VENTRICLE_SURFACE:
00099 case HeartRegionCode::LEFT_VENTRICLE_WALL:
00100 case HeartRegionCode::LEFT_SEPTUM:
00101 lv++;
00102 break;
00103
00104 case HeartRegionCode::RIGHT_VENTRICLE_SURFACE:
00105 case HeartRegionCode::RIGHT_VENTRICLE_WALL:
00106 case HeartRegionCode::RIGHT_SEPTUM:
00107 rv++;
00108 break;
00109
00110 case HeartRegionCode::UNKNOWN:
00111 default:
00112 NEVER_REACHED;
00113 }
00114 }
00115
00116
00117 if (rv>lv)
00118 {
00119 return M_PI/4.0;
00120 }
00121
00122
00123 return M_PI/3.0;
00124 }
00125
00126 template<unsigned SPACE_DIM>
00127 StreeterFibreGenerator<SPACE_DIM>::StreeterFibreGenerator(TetrahedralMesh<SPACE_DIM,SPACE_DIM>& rMesh)
00128 : mrMesh(rMesh),
00129 mpGeometryInfo(NULL)
00130 {
00131 }
00132
00133 template<unsigned SPACE_DIM>
00134 StreeterFibreGenerator<SPACE_DIM>::~StreeterFibreGenerator()
00135 {
00136 delete mpGeometryInfo;
00137 }
00138
00139 template<unsigned SPACE_DIM>
00140 void StreeterFibreGenerator<SPACE_DIM>::SetSurfaceFiles(
00141 const std::string &epicardiumFile,
00142 const std::string &rightVentricleFile,
00143 const std::string &leftVentricleFile)
00144 {
00145
00146 mpGeometryInfo = new HeartGeometryInformation<SPACE_DIM>(mrMesh, epicardiumFile, leftVentricleFile, rightVentricleFile);
00147 }
00148
00149 template<unsigned SPACE_DIM>
00150 void StreeterFibreGenerator<SPACE_DIM>::GenerateOrthotropicFibreOrientation(
00151 std::string outputDirectory,
00152 std::string fibreOrientationFile,
00153 bool logInfo)
00154 {
00155 if (mpGeometryInfo == NULL)
00156 {
00157 EXCEPTION("Files defining the heart surfaces not set");
00158 }
00159
00160
00161 OutputFileHandler handler(outputDirectory, false);
00162 out_stream p_fibre_file = handler.OpenOutputFile(fibreOrientationFile);
00163 out_stream p_regions_file, p_thickness_file, p_ave_thickness_file, p_grad_thickness_file;
00164
00165 if (logInfo)
00166 {
00167 p_regions_file = handler.OpenOutputFile("node_regions.data");
00168 p_thickness_file = handler.OpenOutputFile("wall_thickness.data");
00169 p_ave_thickness_file = handler.OpenOutputFile("averaged_thickness.data");
00170 p_grad_thickness_file = handler.OpenOutputFile("grad_thickness.data");
00171 }
00172
00173
00174 unsigned num_elements = mrMesh.GetNumElements();
00175 *p_fibre_file << num_elements << std::endl;
00176
00177
00178
00179 CheckVentricleAlignment();
00180
00181
00182 unsigned num_nodes = mrMesh.GetNumNodes();
00183 std::vector<double> wall_thickness(num_nodes);
00184 for (unsigned node_index=0; node_index<num_nodes; node_index++)
00185 {
00186 double dist_epi, dist_endo;
00187
00188 HeartRegionType node_region = mpGeometryInfo->GetHeartRegion(node_index);
00189
00190 switch(node_region)
00191 {
00192 case HeartRegionCode::LEFT_VENTRICLE_SURFACE:
00193 case HeartRegionCode::LEFT_VENTRICLE_WALL:
00194 dist_epi = mpGeometryInfo->rGetDistanceMapEpicardium()[node_index];
00195 dist_endo = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
00196 break;
00197
00198 case HeartRegionCode::RIGHT_VENTRICLE_SURFACE:
00199 case HeartRegionCode::RIGHT_VENTRICLE_WALL:
00200 dist_epi = mpGeometryInfo->rGetDistanceMapEpicardium()[node_index];
00201 dist_endo = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
00202 break;
00203
00204 case HeartRegionCode::LEFT_SEPTUM:
00205 dist_epi = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
00206 dist_endo = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
00207 break;
00208
00209 case HeartRegionCode::RIGHT_SEPTUM:
00210 dist_epi = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
00211 dist_endo = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
00212 break;
00213
00214 case HeartRegionCode::UNKNOWN:
00215 #define COVERAGE_IGNORE
00216 std::cerr << "Wrong distances node: " << node_index << "\t"
00217 << "Epi " << mpGeometryInfo->rGetDistanceMapEpicardium()[node_index] << "\t"
00218 << "RV " << mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index] << "\t"
00219 << "LV " << mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index]
00220 << std::endl;
00221
00222
00223 dist_epi = 1;
00224 dist_endo = 0;
00225 break;
00226 #undef COVERAGE_IGNORE
00227
00228 default:
00229 NEVER_REACHED;
00230 }
00231
00232 wall_thickness[node_index] = dist_endo / (dist_endo + dist_epi);
00233
00234 if (std::isnan(wall_thickness[node_index]))
00235 {
00236 #define COVERAGE_IGNORE
00237
00238
00239
00240
00241 wall_thickness[node_index] = 0;
00242 #undef COVERAGE_IGNORE
00243 }
00244
00245 if (logInfo)
00246 {
00247 *p_regions_file << node_region*100 << "\n";
00248 *p_thickness_file << wall_thickness[node_index] << "\n";
00249 }
00250 }
00251
00252
00253
00254
00255 std::vector<double> averaged_wall_thickness(num_nodes);
00256 for (unsigned node_index=0; node_index<num_nodes; node_index++)
00257 {
00258 averaged_wall_thickness[node_index] = GetAveragedThickness(node_index, wall_thickness);
00259
00260 if (logInfo)
00261 {
00262 *p_ave_thickness_file << averaged_wall_thickness[node_index] << "\n";
00263 }
00264
00265 }
00266
00267
00268
00269
00270 c_vector<double,SPACE_DIM> grad_ave_wall_thickness;
00271
00272 for (unsigned element_index=0; element_index<num_elements; element_index++)
00273 {
00274 Element<SPACE_DIM,SPACE_DIM>* p_element = mrMesh.GetElement(element_index);
00275
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291 c_vector<double, SPACE_DIM+1> elem_nodes_ave_thickness;
00292 double element_averaged_thickness = 0.0;
00293 c_vector<HeartRegionType, SPACE_DIM+1> elem_nodes_region;
00294
00295 for (unsigned local_node_index=0; local_node_index<SPACE_DIM+1; local_node_index++)
00296 {
00297
00298 unsigned global_node_index = p_element->GetNode(local_node_index)->GetIndex();
00299
00300 elem_nodes_ave_thickness[local_node_index] = averaged_wall_thickness[global_node_index];
00301 elem_nodes_region[local_node_index] = mpGeometryInfo->GetHeartRegion(global_node_index);
00302
00303
00304 element_averaged_thickness += wall_thickness[global_node_index];
00305 }
00306
00307 element_averaged_thickness /= SPACE_DIM+1;
00308
00310 assert (SPACE_DIM==3);
00311
00312 c_matrix<double, SPACE_DIM+1, SPACE_DIM> basis_functions( zero_matrix<double>(4u,3u) );
00313 basis_functions(0,0) = basis_functions(0,1) = basis_functions(0,2) = -1.0;
00314 basis_functions(1,0) = basis_functions(2,1) = basis_functions(3,2) = 1.0;
00315
00316 c_matrix<double, SPACE_DIM+1, SPACE_DIM> temp;
00317 c_matrix<double, SPACE_DIM, SPACE_DIM> jacobian, inverse_jacobian;
00318 double jacobian_det;
00319 mrMesh.GetInverseJacobianForElement(element_index, jacobian, jacobian_det, inverse_jacobian);
00320 noalias(temp) = prod (basis_functions, inverse_jacobian);
00321 noalias(grad_ave_wall_thickness) = prod(elem_nodes_ave_thickness, temp);
00322
00323 grad_ave_wall_thickness /= norm_2(grad_ave_wall_thickness);
00324
00325 if (logInfo)
00326 {
00327 *p_grad_thickness_file << grad_ave_wall_thickness[0] << " " << grad_ave_wall_thickness[1] << " " << grad_ave_wall_thickness[2] << std::endl;
00328 }
00329
00330
00331
00332
00333
00334
00335
00336
00337 c_vector<double, SPACE_DIM> fibre_direction = VectorProduct(grad_ave_wall_thickness, Create_c_vector(1.0, 0.0, 0.0));
00338 fibre_direction /= norm_2(fibre_direction);
00339
00340
00341
00342
00343 c_vector<double, SPACE_DIM> longitude_direction = VectorProduct(grad_ave_wall_thickness, fibre_direction);
00344
00345
00346
00347
00348
00349
00350
00351
00352 double alpha = GetFibreMaxAngle(elem_nodes_region) * pow( (1 - 2*element_averaged_thickness), 3 );
00353
00354
00355
00356
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378 c_vector<double, SPACE_DIM> rotated_fibre_direction = fibre_direction*cos(alpha) + longitude_direction*sin(alpha);
00379 c_vector<double, SPACE_DIM> rotated_longitude_direction = -fibre_direction*sin(alpha) + longitude_direction*cos(alpha);
00380
00381
00382
00383
00384
00385 assert( fabs(norm_2(rotated_fibre_direction) - 1) < 100*DBL_EPSILON );
00386 assert( fabs(norm_2(grad_ave_wall_thickness) - 1) < 100*DBL_EPSILON );
00387 assert( fabs(norm_2(rotated_longitude_direction) - 1) < 100*DBL_EPSILON );
00388
00389 assert( fabs(inner_prod(rotated_fibre_direction, grad_ave_wall_thickness)) < 100*DBL_EPSILON );
00390 assert( fabs(inner_prod(rotated_fibre_direction, rotated_longitude_direction)) < 100*DBL_EPSILON);
00391 assert( fabs(inner_prod(grad_ave_wall_thickness, rotated_longitude_direction)) < 100*DBL_EPSILON);
00392
00393
00394
00395
00396
00397
00398
00399 *p_fibre_file << rotated_fibre_direction[0] << " " << rotated_fibre_direction[1] << " " << rotated_fibre_direction[2] << " "
00400 << grad_ave_wall_thickness[0] << " " << grad_ave_wall_thickness[1] << " " << grad_ave_wall_thickness[2] << " "
00401 << rotated_longitude_direction[0] << " " << rotated_longitude_direction[1] << " " << rotated_longitude_direction[2] << std::endl;
00402 }
00403
00404 p_fibre_file->close();
00405
00406 if (logInfo)
00407 {
00408 p_regions_file->close();
00409 p_thickness_file->close();
00410 p_ave_thickness_file->close();
00411 p_grad_thickness_file->close();
00412 }
00413 }
00414
00415 template<unsigned SPACE_DIM>
00416 void StreeterFibreGenerator<SPACE_DIM>::CheckVentricleAlignment()
00417 {
00418 double min_y_rv=DBL_MAX;
00419 double max_y_rv=-DBL_MAX;
00420 double average_y_rv=0.0;
00421 for (unsigned i=0; i<mpGeometryInfo->rGetNodesOnRVSurface().size(); i++)
00422 {
00423 double y=mrMesh.GetNode(mpGeometryInfo->rGetNodesOnRVSurface()[i])->rGetLocation()[1];
00424 average_y_rv += y;
00425 if (y<min_y_rv)
00426 {
00427 min_y_rv=y;
00428 }
00429
00430 if (y>max_y_rv)
00431 {
00432 max_y_rv=y;
00433 }
00434 }
00435 average_y_rv /= mpGeometryInfo->rGetNodesOnRVSurface().size();
00436
00437 double min_y_lv=DBL_MAX;
00438 double max_y_lv=-DBL_MAX;
00439 double average_y_lv=0.0;
00440 for (unsigned i=0; i<mpGeometryInfo->rGetNodesOnLVSurface().size(); i++)
00441 {
00442 double y=mrMesh.GetNode(mpGeometryInfo->rGetNodesOnLVSurface()[i])->rGetLocation()[1];
00443 average_y_lv += y;
00444 if (y<min_y_lv)
00445 {
00446 min_y_lv=y;
00447 }
00448
00449 if (y>max_y_lv)
00450 {
00451 max_y_lv=y;
00452 }
00453 }
00454 average_y_lv /= mpGeometryInfo->rGetNodesOnLVSurface().size();
00455
00456
00457
00458
00459
00460 assert(average_y_lv < min_y_rv || average_y_lv > max_y_rv);
00461
00462
00463 assert(average_y_rv < min_y_lv || average_y_rv > max_y_lv);
00464 }
00465
00467
00469
00470
00471
00472 template class StreeterFibreGenerator<3>;