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 "StreeterFibreGenerator.hpp"
00030
00031 #include <fstream>
00032 #include <sstream>
00033 #include "OutputFileHandler.hpp"
00034 #include "Exception.hpp"
00035
00036
00037 template<unsigned SPACE_DIM>
00038 typename StreeterFibreGenerator<SPACE_DIM>::RegionType_
00039 StreeterFibreGenerator<SPACE_DIM>::GetHeartRegion(unsigned nodeIndex) const
00040 {
00041
00042 if (mDistMapRightVentricle[nodeIndex] >= mDistMapEpicardium[nodeIndex] &&
00043 mDistMapRightVentricle[nodeIndex] >= mDistMapLeftVentricle[nodeIndex])
00044 {
00045 return LEFT_VENTRICLE_WALL;
00046 }
00047
00048 if (mDistMapLeftVentricle[nodeIndex] >= mDistMapEpicardium[nodeIndex] &&
00049 mDistMapLeftVentricle[nodeIndex] >= mDistMapRightVentricle[nodeIndex])
00050 {
00051 return RIGHT_VENTRICLE_WALL;
00052 }
00053
00054 if (mDistMapEpicardium[nodeIndex] >= mDistMapLeftVentricle[nodeIndex] &&
00055 mDistMapEpicardium[nodeIndex] >= mDistMapRightVentricle[nodeIndex])
00056 {
00057 if (mDistMapLeftVentricle[nodeIndex]
00058 < LEFT_SEPTUM_SIZE*(mDistMapLeftVentricle[nodeIndex] + mDistMapRightVentricle[nodeIndex]))
00059 {
00060 return LEFT_SEPTUM;
00061 }
00062 else
00063 {
00064 return RIGHT_SEPTUM;
00065 }
00066 }
00067
00068 return UNKNOWN;
00069 }
00070
00071
00072 template<unsigned SPACE_DIM>
00073 double StreeterFibreGenerator<SPACE_DIM>::GetAveragedThickness(
00074 const unsigned nodeIndex, const std::vector<double>& wallThickness) const
00075 {
00076
00077 double average = wallThickness[nodeIndex];
00078 unsigned nodes_visited = 1;
00079
00080
00081 std::set<unsigned> visited_nodes;
00082 visited_nodes.insert(nodeIndex);
00083
00084 Node<SPACE_DIM>* p_current_node = mrMesh.GetNode(nodeIndex);
00085
00087
00088
00089 for(typename Node<SPACE_DIM>::ContainingElementIterator element_iterator = p_current_node->ContainingElementsBegin();
00090 element_iterator != p_current_node->ContainingElementsEnd();
00091 ++element_iterator)
00092 {
00093
00094 Element<SPACE_DIM,SPACE_DIM>* p_containing_element = mrMesh.GetElement(*element_iterator);
00095
00096
00097 for(unsigned node_local_index=0;
00098 node_local_index<p_containing_element->GetNumNodes();
00099 node_local_index++)
00100 {
00101 Node<SPACE_DIM>* p_neighbour_node = p_containing_element->GetNode(node_local_index);
00102 unsigned neighbour_node_index = p_neighbour_node->GetIndex();
00103
00104
00105 if (visited_nodes.find(neighbour_node_index) == visited_nodes.end())
00106 {
00107 average += wallThickness[neighbour_node_index];
00108 visited_nodes.insert(neighbour_node_index);
00109 nodes_visited++;
00110 }
00111 }
00112 }
00113
00114 return average/nodes_visited;
00115 }
00116
00117
00118 template<unsigned SPACE_DIM>
00119 void StreeterFibreGenerator<SPACE_DIM>::ProcessLine(
00120 const std::string& line, std::set<unsigned>& surfaceNodeIndexSet) const
00121 {
00122 unsigned num_nodes = 0;
00123 std::stringstream line_stream(line);
00124
00125 while (!line_stream.eof())
00126 {
00127 unsigned item;
00128 line_stream >> item;
00129
00130 surfaceNodeIndexSet.insert(item-1);
00131
00132 num_nodes++;
00133 }
00134
00135 if (num_nodes != SPACE_DIM)
00136 {
00137 EXCEPTION("Wrong file format");
00138 }
00139
00140 }
00141
00142
00143 template<unsigned SPACE_DIM>
00144 void StreeterFibreGenerator<SPACE_DIM>::GetNodesAtSurface(
00145 const std::string& surfaceFile, std::vector<unsigned>& surfaceVector) const
00146 {
00147
00148 std::ifstream file_stream;
00149 file_stream.open(surfaceFile.c_str());
00150 if (!file_stream.is_open())
00151 {
00152 EXCEPTION("Wrong surface definition file name.");
00153 }
00154
00155
00156 std::set<unsigned> surface_node_index_set;
00157
00158
00159 std::string line;
00160 getline(file_stream, line);
00161 do
00162 {
00163 ProcessLine(line, surface_node_index_set);
00164
00165 getline(file_stream, line);
00166 }
00167 while(!file_stream.eof());
00168
00169
00170 surfaceVector.reserve(surface_node_index_set.size());
00171
00172
00173 for(std::set<unsigned>::iterator node_index_it=surface_node_index_set.begin();
00174 node_index_it != surface_node_index_set.end();
00175 node_index_it++)
00176 {
00177 surfaceVector.push_back(*node_index_it);
00178 }
00179
00180 file_stream.close();
00181 }
00182
00183
00184 template<unsigned SPACE_DIM>
00185 double StreeterFibreGenerator<SPACE_DIM>::GetFibreMaxAngle(
00186 const c_vector<RegionType_, SPACE_DIM+1>& nodesRegion) const
00187 {
00188 unsigned lv=0, rv=0;
00189
00190 for (unsigned index=0; index<SPACE_DIM+1; index++)
00191 {
00192 switch (nodesRegion[index])
00193 {
00194 case LEFT_VENTRICLE_WALL:
00195 case LEFT_SEPTUM:
00196 lv++;
00197 break;
00198
00199 case RIGHT_VENTRICLE_WALL:
00200 case RIGHT_SEPTUM:
00201 rv++;
00202 break;
00203
00204 case UNKNOWN:
00205 break;
00206 }
00207 }
00208
00209
00210 if (rv>lv)
00211 {
00212 return M_PI/4.0;
00213 }
00214
00215
00216 return M_PI/3.0;
00217 }
00218
00219 template<unsigned SPACE_DIM>
00220 StreeterFibreGenerator<SPACE_DIM>::StreeterFibreGenerator(TetrahedralMesh<SPACE_DIM,SPACE_DIM>& rMesh)
00221 : mrMesh(rMesh),
00222 mFilesSet(false)
00223 {
00224 mNumNodes = mrMesh.GetNumNodes();
00225 mNumElements = mrMesh.GetNumElements();
00226 mpDistanceCalculator = new DistanceMapCalculator<SPACE_DIM>(mrMesh);
00227 }
00228
00229 template<unsigned SPACE_DIM>
00230 StreeterFibreGenerator<SPACE_DIM>::~StreeterFibreGenerator()
00231 {
00232 delete mpDistanceCalculator;
00233 }
00234
00235 template<unsigned SPACE_DIM>
00236 void StreeterFibreGenerator<SPACE_DIM>::SetSurfaceFiles(
00237 std::string epicardiumFile,
00238 std::string rightVentricleFile,
00239 std::string leftVentricleFile)
00240 {
00241 mEpiFile = epicardiumFile;
00242 mRVFile = rightVentricleFile;
00243 mLVFile = leftVentricleFile;
00244
00245 mFilesSet = true;
00246 }
00247
00248 template<unsigned SPACE_DIM>
00249 void StreeterFibreGenerator<SPACE_DIM>::GenerateOrthotropicFibreOrientation(
00250 std::string outputDirectory,
00251 std::string fibreOrientationFile,
00252 bool logInfo)
00253 {
00254 if (!mFilesSet)
00255 {
00256 EXCEPTION("Files defining the heart surfaces not set");
00257 }
00258
00259
00260 OutputFileHandler handler(outputDirectory, false);
00261 out_stream p_fibre_file = handler.OpenOutputFile(fibreOrientationFile);
00262 out_stream p_regions_file, p_thickness_file, p_ave_thickness_file, p_grad_thickness_file;
00263
00264 if (logInfo)
00265 {
00266 p_regions_file = handler.OpenOutputFile("node_regions.data");
00267 p_thickness_file = handler.OpenOutputFile("wall_thickness.data");
00268 p_ave_thickness_file = handler.OpenOutputFile("averaged_thickness.data");
00269 p_grad_thickness_file = handler.OpenOutputFile("grad_thickness.data");
00270 }
00271
00272
00273 *p_fibre_file << mNumElements << std::endl;
00274
00275
00276 GetNodesAtSurface(mEpiFile, mEpiSurface);
00277 GetNodesAtSurface(mRVFile, mRVSurface);
00278 GetNodesAtSurface(mLVFile, mLVSurface);
00279
00280 CheckVentricleAlignment();
00281
00282
00283 mpDistanceCalculator->ComputeDistanceMap(mEpiSurface, mDistMapEpicardium);
00284 mpDistanceCalculator->ComputeDistanceMap(mRVSurface, mDistMapRightVentricle);
00285 mpDistanceCalculator->ComputeDistanceMap(mLVSurface, mDistMapLeftVentricle);
00286
00287
00288
00289
00290 std::vector<double> wall_thickness(mNumNodes);
00291 for (unsigned node_index=0; node_index<mNumNodes; node_index++)
00292 {
00293 double dist_epi, dist_endo;
00294
00295 RegionType_ node_region = GetHeartRegion(node_index);
00296
00297 switch(node_region)
00298 {
00299 case LEFT_VENTRICLE_WALL:
00300 dist_epi = mDistMapEpicardium[node_index];
00301 dist_endo = mDistMapLeftVentricle[node_index];
00302 break;
00303
00304 case RIGHT_VENTRICLE_WALL:
00305 dist_epi = mDistMapEpicardium[node_index];
00306 dist_endo = mDistMapRightVentricle[node_index];
00307 break;
00308
00309 case LEFT_SEPTUM:
00310 dist_epi = mDistMapRightVentricle[node_index];
00311 dist_endo = mDistMapLeftVentricle[node_index];
00312 break;
00313
00314 case RIGHT_SEPTUM:
00315 dist_epi = mDistMapLeftVentricle[node_index];
00316 dist_endo = mDistMapRightVentricle[node_index];
00317 break;
00318
00319 case UNKNOWN:
00320 #define COVERAGE_IGNORE
00321 std::cerr << "Wrong distances node: " << node_index << "\t"
00322 << "Epi " << mDistMapEpicardium[node_index] << "\t"
00323 << "RV " << mDistMapRightVentricle[node_index] << "\t"
00324 << "LV " << mDistMapLeftVentricle[node_index]
00325 << std::endl;
00326
00327
00328 dist_epi = 1;
00329 dist_endo = 0;
00330 break;
00331 #undef COVERAGE_IGNORE
00332 }
00333
00334 wall_thickness[node_index] = dist_endo / (dist_endo + dist_epi);
00335
00336 if (isnan(wall_thickness[node_index]))
00337 {
00338 #define COVERAGE_IGNORE
00339
00340
00341
00342
00343 wall_thickness[node_index] = 0;
00344 #undef COVERAGE_IGNORE
00345 }
00346
00347 if (logInfo)
00348 {
00349 *p_regions_file << node_region*100 << "\n";
00350 *p_thickness_file << wall_thickness[node_index] << "\n";
00351 }
00352 }
00353
00354
00355
00356
00357 std::vector<double> averaged_wall_thickness(mNumNodes);
00358 for (unsigned node_index=0; node_index<mNumNodes; node_index++)
00359 {
00360 averaged_wall_thickness[node_index] = GetAveragedThickness(node_index, wall_thickness);
00361
00362 if (logInfo)
00363 {
00364 *p_ave_thickness_file << averaged_wall_thickness[node_index] << "\n";
00365 }
00366
00367 }
00368
00369
00370
00371
00372 c_vector<double,SPACE_DIM> grad_ave_wall_thickness;
00373
00374 for (unsigned element_index=0; element_index<mNumElements; element_index++)
00375 {
00376 Element<SPACE_DIM,SPACE_DIM>* p_element = mrMesh.GetElement(element_index);
00377
00378
00379
00380
00381
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393 c_vector<double, SPACE_DIM+1> elem_nodes_ave_thickness;
00394 double element_averaged_thickness = 0.0;
00395 c_vector<RegionType_, SPACE_DIM+1> elem_nodes_region;
00396
00397 for (unsigned local_node_index=0; local_node_index<SPACE_DIM+1; local_node_index++)
00398 {
00399
00400 unsigned global_node_index = p_element->GetNode(local_node_index)->GetIndex();
00401
00402 elem_nodes_ave_thickness[local_node_index] = averaged_wall_thickness[global_node_index];
00403 elem_nodes_region[local_node_index] = GetHeartRegion(global_node_index);
00404
00405
00406 element_averaged_thickness += wall_thickness[global_node_index];
00407 }
00408
00409 element_averaged_thickness /= SPACE_DIM+1;
00410
00412 assert (SPACE_DIM==3);
00413
00414 c_matrix<double, SPACE_DIM+1, SPACE_DIM> basis_functions( zero_matrix<double>(4u,3u) );
00415 basis_functions(0,0) = basis_functions(0,1) = basis_functions(0,2) = -1.0;
00416 basis_functions(1,0) = basis_functions(2,1) = basis_functions(3,2) = 1.0;
00417
00418 c_matrix<double, SPACE_DIM+1, SPACE_DIM> temp;
00419 c_matrix<double, SPACE_DIM, SPACE_DIM> jacobian, inverse_jacobian;
00420 double jacobian_det;
00421 mrMesh.GetInverseJacobianForElement(element_index, jacobian, jacobian_det, inverse_jacobian);
00422 noalias(temp) = prod (basis_functions, inverse_jacobian);
00423 noalias(grad_ave_wall_thickness) = prod(elem_nodes_ave_thickness, temp);
00424
00425 grad_ave_wall_thickness /= norm_2(grad_ave_wall_thickness);
00426
00427 if (logInfo)
00428 {
00429 *p_grad_thickness_file << grad_ave_wall_thickness[0] << " " << grad_ave_wall_thickness[1] << " " << grad_ave_wall_thickness[2] << std::endl;
00430 }
00431
00432
00433
00434
00435
00436
00437
00438
00439 c_vector<double, SPACE_DIM> fibre_direction = VectorProduct(grad_ave_wall_thickness, Create_c_vector(1.0, 0.0, 0.0));
00440 fibre_direction /= norm_2(fibre_direction);
00441
00442
00443
00444
00445 c_vector<double, SPACE_DIM> longitude_direction = VectorProduct(grad_ave_wall_thickness, fibre_direction);
00446
00447
00448
00449
00450
00451
00452
00453
00454 double alpha = GetFibreMaxAngle(elem_nodes_region) * pow( (1 - 2*element_averaged_thickness), 3 );
00455
00456
00457
00458
00459
00460
00461
00462
00463
00464
00465
00466
00467
00468
00469
00470
00471
00472
00473
00474
00475
00476
00477
00478
00479
00480 c_vector<double, SPACE_DIM> rotated_fibre_direction = fibre_direction*cos(alpha) + longitude_direction*sin(alpha);
00481 c_vector<double, SPACE_DIM> rotated_longitude_direction = -fibre_direction*sin(alpha) + longitude_direction*cos(alpha);
00482
00483
00484
00485
00486
00487 assert( fabs(norm_2(rotated_fibre_direction) - 1) < 100*DBL_EPSILON );
00488 assert( fabs(norm_2(grad_ave_wall_thickness) - 1) < 100*DBL_EPSILON );
00489 assert( fabs(norm_2(rotated_longitude_direction) - 1) < 100*DBL_EPSILON );
00490
00491 assert( fabs(inner_prod(rotated_fibre_direction, grad_ave_wall_thickness)) < 100*DBL_EPSILON );
00492 assert( fabs(inner_prod(rotated_fibre_direction, rotated_longitude_direction)) < 100*DBL_EPSILON);
00493 assert( fabs(inner_prod(grad_ave_wall_thickness, rotated_longitude_direction)) < 100*DBL_EPSILON);
00494
00495
00496
00497
00498
00499
00500
00501 *p_fibre_file << rotated_fibre_direction[0] << " " << rotated_fibre_direction[1] << " " << rotated_fibre_direction[2] << " "
00502 << grad_ave_wall_thickness[0] << " " << grad_ave_wall_thickness[1] << " " << grad_ave_wall_thickness[2] << " "
00503 << rotated_longitude_direction[0] << " " << rotated_longitude_direction[1] << " " << rotated_longitude_direction[2] << std::endl;
00504 }
00505
00506 p_fibre_file->close();
00507
00508 if (logInfo)
00509 {
00510 p_regions_file->close();
00511 p_thickness_file->close();
00512 p_ave_thickness_file->close();
00513 p_grad_thickness_file->close();
00514 }
00515 }
00516
00517 template<unsigned SPACE_DIM>
00518 void StreeterFibreGenerator<SPACE_DIM>::CheckVentricleAlignment()
00519 {
00520 double min_y_rv=DBL_MAX;
00521 double max_y_rv=-DBL_MAX;
00522 double average_y_rv=0.0;
00523 for (unsigned i=0; i<mRVSurface.size(); i++)
00524 {
00525 double y=mrMesh.GetNode(mRVSurface[i])->rGetLocation()[1];
00526 average_y_rv += y;
00527 if (y<min_y_rv)
00528 {
00529 min_y_rv=y;
00530 }
00531
00532 if (y>max_y_rv)
00533 {
00534 max_y_rv=y;
00535 }
00536 }
00537 average_y_rv /= mRVSurface.size();
00538
00539 double min_y_lv=DBL_MAX;
00540 double max_y_lv=-DBL_MAX;
00541 double average_y_lv=0.0;
00542 for (unsigned i=0; i<mLVSurface.size(); i++)
00543 {
00544 double y=mrMesh.GetNode(mLVSurface[i])->rGetLocation()[1];
00545 average_y_lv += y;
00546 if (y<min_y_lv)
00547 {
00548 min_y_lv=y;
00549 }
00550
00551 if (y>max_y_lv)
00552 {
00553 max_y_lv=y;
00554 }
00555 }
00556 average_y_lv /= mLVSurface.size();
00557
00558
00559
00560
00561
00562 assert(average_y_lv < min_y_rv || average_y_lv > max_y_rv);
00563
00564
00565 assert(average_y_rv < min_y_lv || average_y_rv > max_y_lv);
00566 }
00567
00569
00571
00572
00573
00574 template class StreeterFibreGenerator<3>;