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
00039
00040 template<unsigned SPACE_DIM>
00041 double StreeterFibreGenerator<SPACE_DIM>::GetAveragedThicknessLocalNode(
00042 const unsigned nodeIndex, const std::vector<double>& wallThickness) const
00043 {
00044 if (nodeIndex < mrMesh.GetDistributedVectorFactory()->GetLow() ||
00045 nodeIndex >= mrMesh.GetDistributedVectorFactory()->GetHigh() )
00046 {
00047 return 0.0;
00048 }
00049
00050
00051 double average = wallThickness[nodeIndex];
00052 unsigned nodes_visited = 1;
00053
00054
00055 std::set<unsigned> visited_nodes;
00056 visited_nodes.insert(nodeIndex);
00057
00058 Node<SPACE_DIM>* p_current_node = mrMesh.GetNode(nodeIndex);
00059
00060
00061 for(typename Node<SPACE_DIM>::ContainingElementIterator element_iterator = p_current_node->ContainingElementsBegin();
00062 element_iterator != p_current_node->ContainingElementsEnd();
00063 ++element_iterator)
00064 {
00065
00066 Element<SPACE_DIM,SPACE_DIM>* p_containing_element = mrMesh.GetElement(*element_iterator);
00067
00068
00069 for(unsigned node_local_index=0;
00070 node_local_index<p_containing_element->GetNumNodes();
00071 node_local_index++)
00072 {
00073 Node<SPACE_DIM>* p_neighbour_node = p_containing_element->GetNode(node_local_index);
00074 unsigned neighbour_node_index = p_neighbour_node->GetIndex();
00075
00076
00077 if (visited_nodes.find(neighbour_node_index) == visited_nodes.end())
00078 {
00079 average += wallThickness[neighbour_node_index];
00080 visited_nodes.insert(neighbour_node_index);
00081 nodes_visited++;
00082 }
00083 }
00084 }
00085
00086 return average/nodes_visited;
00087 }
00088
00089
00090
00091
00092 template<unsigned SPACE_DIM>
00093 double StreeterFibreGenerator<SPACE_DIM>::GetFibreMaxAngle(
00094 const c_vector<HeartRegionType, SPACE_DIM+1>& nodesRegionsForElement) const
00095 {
00096 unsigned lv=0, rv=0;
00097
00098 for (unsigned index=0; index<SPACE_DIM+1; index++)
00099 {
00100 switch (nodesRegionsForElement[index])
00101 {
00102 case HeartGeometryInformation<SPACE_DIM>::LEFT_VENTRICLE_SURFACE:
00103 case HeartGeometryInformation<SPACE_DIM>::LEFT_VENTRICLE_WALL:
00104 case HeartGeometryInformation<SPACE_DIM>::LEFT_SEPTUM:
00105 lv++;
00106 break;
00107
00108 case HeartGeometryInformation<SPACE_DIM>::RIGHT_VENTRICLE_SURFACE:
00109 case HeartGeometryInformation<SPACE_DIM>::RIGHT_VENTRICLE_WALL:
00110 case HeartGeometryInformation<SPACE_DIM>::RIGHT_SEPTUM:
00111 rv++;
00112 break;
00113
00114 case HeartGeometryInformation<SPACE_DIM>::UNKNOWN:
00115 default:
00116 NEVER_REACHED;
00117 }
00118 }
00119
00120
00121 if (rv>lv)
00122 {
00123 return M_PI/4.0;
00124 }
00125
00126
00127 return M_PI/3.0;
00128 }
00129
00130 template<unsigned SPACE_DIM>
00131 StreeterFibreGenerator<SPACE_DIM>::StreeterFibreGenerator(AbstractTetrahedralMesh<SPACE_DIM,SPACE_DIM>& rMesh)
00132 : mrMesh(rMesh),
00133 mpGeometryInfo(NULL),
00134 mApexToBase(zero_vector<double>(SPACE_DIM))
00135 {
00136 }
00137
00138 template<unsigned SPACE_DIM>
00139 StreeterFibreGenerator<SPACE_DIM>::~StreeterFibreGenerator()
00140 {
00141 delete mpGeometryInfo;
00142 }
00143
00144 template<unsigned SPACE_DIM>
00145 void StreeterFibreGenerator<SPACE_DIM>::SetSurfaceFiles(
00146 const std::string &rEpicardiumFile,
00147 const std::string &rRightVentricleFile,
00148 const std::string &rLeftVentricleFile,
00149 bool indexFromZero)
00150 {
00151
00152 mpGeometryInfo = new HeartGeometryInformation<SPACE_DIM>(mrMesh, rEpicardiumFile, rLeftVentricleFile, rRightVentricleFile, indexFromZero);
00153 }
00154
00155 template<unsigned SPACE_DIM>
00156 void StreeterFibreGenerator<SPACE_DIM>::GenerateOrthotropicFibreOrientation(
00157 std::string outputDirectory,
00158 std::string fibreOrientationFile,
00159 bool logInfo)
00160 {
00161 assert(SPACE_DIM == 3);
00162 if (mpGeometryInfo == NULL)
00163 {
00164 EXCEPTION("Files defining the heart surfaces not set");
00165 }
00166 unsigned num_elements = mrMesh.GetNumElements();
00167
00168
00169 OutputFileHandler handler(outputDirectory, false);
00170 out_stream p_fibre_file, p_regions_file, p_thickness_file, p_ave_thickness_file;
00171
00172
00173
00174 if (PetscTools::AmMaster())
00175 {
00176
00177 p_fibre_file=handler.OpenOutputFile(fibreOrientationFile);
00178
00179 *p_fibre_file << num_elements << std::endl;
00180 p_fibre_file->close();
00181 }
00182 else
00183 {
00184 logInfo=false;
00185 }
00186
00187 if (logInfo)
00188 {
00189 p_regions_file = handler.OpenOutputFile("node_regions.data");
00190 p_thickness_file = handler.OpenOutputFile("wall_thickness.data");
00191 p_ave_thickness_file = handler.OpenOutputFile("averaged_thickness.data");
00192 }
00193
00194
00195
00196 if (fabs(norm_2(mApexToBase)) < DBL_EPSILON)
00197 {
00198 EXCEPTION("Apex to base vector has not been set");
00199 }
00200
00201
00202 unsigned num_nodes = mrMesh.GetNumNodes();
00203 std::vector<double> wall_thickness(num_nodes);
00204 for (unsigned node_index=0; node_index<num_nodes; node_index++)
00205 {
00206 double dist_epi, dist_endo;
00207
00208 HeartRegionType node_region = mpGeometryInfo->GetHeartRegion(node_index);
00209
00210 switch(node_region)
00211 {
00212 case HeartGeometryInformation<SPACE_DIM>::LEFT_VENTRICLE_SURFACE:
00213 case HeartGeometryInformation<SPACE_DIM>::LEFT_VENTRICLE_WALL:
00214 dist_epi = mpGeometryInfo->rGetDistanceMapEpicardium()[node_index];
00215 dist_endo = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
00216 break;
00217
00218 case HeartGeometryInformation<SPACE_DIM>::RIGHT_VENTRICLE_SURFACE:
00219 case HeartGeometryInformation<SPACE_DIM>::RIGHT_VENTRICLE_WALL:
00220 dist_epi = mpGeometryInfo->rGetDistanceMapEpicardium()[node_index];
00221 dist_endo = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
00222 break;
00223
00224 case HeartGeometryInformation<SPACE_DIM>::LEFT_SEPTUM:
00225 dist_epi = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
00226 dist_endo = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
00227 break;
00228
00229 case HeartGeometryInformation<SPACE_DIM>::RIGHT_SEPTUM:
00230 dist_epi = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
00231 dist_endo = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
00232 break;
00233
00234 case HeartGeometryInformation<SPACE_DIM>::UNKNOWN:
00235 #define COVERAGE_IGNORE
00236 std::cerr << "Wrong distances node: " << node_index << "\t"
00237 << "Epi " << mpGeometryInfo->rGetDistanceMapEpicardium()[node_index] << "\t"
00238 << "RV " << mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index] << "\t"
00239 << "LV " << mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index]
00240 << std::endl;
00241
00242
00243 dist_epi = 1;
00244 dist_endo = 0;
00245 break;
00246 #undef COVERAGE_IGNORE
00247
00248 default:
00249 NEVER_REACHED;
00250 }
00251
00252 wall_thickness[node_index] = dist_endo / (dist_endo + dist_epi);
00253
00254 if (std::isnan(wall_thickness[node_index]))
00255 {
00256 #define COVERAGE_IGNORE
00257
00258
00259
00260
00261 wall_thickness[node_index] = 0;
00262 #undef COVERAGE_IGNORE
00263 }
00264
00265 if (logInfo)
00266 {
00267 *p_regions_file << node_region*100 << "\n";
00268 *p_thickness_file << wall_thickness[node_index] << "\n";
00269 }
00270 }
00271
00272
00273
00274
00275 std::vector<double> my_averaged_wall_thickness(num_nodes);
00276 std::vector<double> averaged_wall_thickness(num_nodes);
00277 for (unsigned node_index=0; node_index<num_nodes; node_index++)
00278 {
00279 my_averaged_wall_thickness[node_index] = GetAveragedThicknessLocalNode(node_index, wall_thickness);
00280 }
00281
00282
00283 MPI_Allreduce(&my_averaged_wall_thickness[0], &averaged_wall_thickness[0], num_nodes,
00284 MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
00285 if (logInfo)
00286 {
00287 for (unsigned node_index=0; node_index<num_nodes; node_index++)
00288 {
00289 *p_ave_thickness_file << averaged_wall_thickness[node_index] << "\n";
00290 }
00291
00292 }
00293
00294
00295
00296
00297 c_vector<double,SPACE_DIM> grad_ave_wall_thickness;
00298
00299 bool fibre_file_is_open= false;
00300 for (unsigned element_index=0; element_index<num_elements; element_index++)
00301 {
00302
00303 if (mrMesh.CalculateDesignatedOwnershipOfElement(element_index))
00304 {
00305
00306
00307
00308
00309
00310
00311 PetscTools::Barrier("GenerateOrthotropicFibreOrientation");
00312 if (fibre_file_is_open != true)
00313 {
00314
00315 p_fibre_file = handler.OpenOutputFile(fibreOrientationFile, std::ios::app);
00316 fibre_file_is_open=true;
00317 }
00318
00319
00320 Element<SPACE_DIM,SPACE_DIM>* p_element = mrMesh.GetElement(element_index);
00321
00322
00323
00324
00325
00326
00327
00328
00329
00330
00331
00332
00333
00334
00335
00336
00337 c_vector<double, SPACE_DIM+1> elem_nodes_ave_thickness;
00338 double element_averaged_thickness = 0.0;
00339 c_vector<HeartRegionType, SPACE_DIM+1> elem_nodes_region;
00340
00341 for (unsigned local_node_index=0; local_node_index<SPACE_DIM+1; local_node_index++)
00342 {
00343
00344 unsigned global_node_index = p_element->GetNode(local_node_index)->GetIndex();
00345
00346 elem_nodes_ave_thickness[local_node_index] = averaged_wall_thickness[global_node_index];
00347 elem_nodes_region[local_node_index] = mpGeometryInfo->GetHeartRegion(global_node_index);
00348
00349
00350 element_averaged_thickness += wall_thickness[global_node_index];
00351 }
00352
00353 element_averaged_thickness /= SPACE_DIM+1;
00354
00355 c_matrix<double, SPACE_DIM+1, SPACE_DIM> basis_functions( zero_matrix<double>(4u,3u) );
00356 basis_functions(0,0) = basis_functions(0,1) = basis_functions(0,2) = -1.0;
00357 basis_functions(1,0) = basis_functions(2,1) = basis_functions(3,2) = 1.0;
00358
00359 c_matrix<double, SPACE_DIM+1, SPACE_DIM> temp;
00360 c_matrix<double, SPACE_DIM, SPACE_DIM> jacobian, inverse_jacobian;
00361 double jacobian_det;
00362 mrMesh.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
00370
00371
00372
00373
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
00380
00381 c_vector<double, SPACE_DIM> longitude_direction = VectorProduct(grad_ave_wall_thickness, fibre_direction);
00382
00383
00384
00385
00386
00387
00388
00389
00390 double alpha = GetFibreMaxAngle(elem_nodes_region) * SmallPow( (1 - 2*element_averaged_thickness), 3 );
00391
00392
00393
00394
00395
00396
00397
00398
00399
00400
00401
00402
00403
00404
00405
00406
00407
00408
00409
00410
00411
00412
00413
00414
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
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
00433
00434
00435
00436
00437 *p_fibre_file << rotated_fibre_direction[0] << " " << rotated_fibre_direction[1] << " " << rotated_fibre_direction[2] << " "
00438 << grad_ave_wall_thickness[0] << " " << grad_ave_wall_thickness[1] << " " << grad_ave_wall_thickness[2] << " "
00439 << rotated_longitude_direction[0] << " " << rotated_longitude_direction[1] << " " << rotated_longitude_direction[2] << std::endl;
00440 }
00441 else
00442 {
00443
00444
00445 if (fibre_file_is_open)
00446 {
00447 p_fibre_file->close();
00448 fibre_file_is_open=false;
00449 }
00450 PetscTools::Barrier("GenerateOrthotropicFibreOrientation");
00451
00452 }
00453 }
00454
00455 if (fibre_file_is_open)
00456 {
00457 p_fibre_file->close();
00458 }
00459
00460 if (logInfo)
00461 {
00462 p_regions_file->close();
00463 p_thickness_file->close();
00464 p_ave_thickness_file->close();
00465 }
00466
00467 }
00468
00469
00470 template<unsigned SPACE_DIM>
00471 void StreeterFibreGenerator<SPACE_DIM>::SetApexToBase(const c_vector<double, SPACE_DIM>& apexToBase)
00472 {
00473 double norm = norm_2(apexToBase);
00474 if (norm < DBL_EPSILON)
00475 {
00476 EXCEPTION("Apex to base vector should be non-zero");
00477 }
00478 mApexToBase = apexToBase / norm;
00479 }
00480 template<unsigned SPACE_DIM>
00481 void StreeterFibreGenerator<SPACE_DIM>::SetApexToBase(unsigned axis)
00482 {
00483 if (axis >= SPACE_DIM)
00484 {
00485 EXCEPTION("Apex to base coordinate axis was out of range");
00486 }
00487 mApexToBase = zero_vector<double>(SPACE_DIM);
00488 mApexToBase[axis] = 1.0;
00489 }
00490
00492
00494 #define COVERAGE_IGNORE
00495 template class StreeterFibreGenerator<3>;
00496 #undef COVERAGE_IGNORE