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
00030
00031
00032
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
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 < this->mpMesh->GetDistributedVectorFactory()->GetLow() ||
00052 nodeIndex >= this->mpMesh->GetDistributedVectorFactory()->GetHigh() )
00053 {
00054 return 0.0;
00055 }
00056
00057
00058 double average = wallThickness[nodeIndex];
00059 unsigned nodes_visited = 1;
00060
00061
00062 std::set<unsigned> visited_nodes;
00063 visited_nodes.insert(nodeIndex);
00064
00065 Node<SPACE_DIM>* p_current_node = this->mpMesh->GetNode(nodeIndex);
00066
00067
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
00073 Element<SPACE_DIM,SPACE_DIM>* p_containing_element = this->mpMesh->GetElement(*element_iterator);
00074
00075
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
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
00128 if (rv>lv)
00129 {
00130 return M_PI/4.0;
00131 }
00132
00133
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 : AbstractPerElementWriter<SPACE_DIM,SPACE_DIM,SPACE_DIM*SPACE_DIM>(&rMesh),
00140 mpGeometryInfo(NULL),
00141 mApexToBase(zero_vector<double>(SPACE_DIM)),
00142 mLogInfo(false)
00143 {
00144 mWallThickness.resize(rMesh.GetNumNodes());
00145 mAveragedWallThickness.resize(rMesh.GetNumNodes());
00146 }
00147
00148 template<unsigned SPACE_DIM>
00149 StreeterFibreGenerator<SPACE_DIM>::~StreeterFibreGenerator()
00150 {
00151 delete mpGeometryInfo;
00152 }
00153
00154 template<unsigned SPACE_DIM>
00155 void StreeterFibreGenerator<SPACE_DIM>::SetSurfaceFiles(
00156 const std::string &rEpicardiumFile,
00157 const std::string &rRightVentricleFile,
00158 const std::string &rLeftVentricleFile,
00159 bool indexFromZero)
00160 {
00161
00162 mpGeometryInfo = new HeartGeometryInformation<SPACE_DIM>(*(this->mpMesh), rEpicardiumFile, rLeftVentricleFile, rRightVentricleFile, indexFromZero);
00163 }
00164
00165 template<unsigned SPACE_DIM>
00166 void StreeterFibreGenerator<SPACE_DIM>::WriteHeaderOnMaster()
00167 {
00168 *(this->mpMasterFile) << this->mpMesh->GetNumElements();
00169 *(this->mpMasterFile) << std::setprecision(16);
00170 }
00171
00172 template<unsigned SPACE_DIM>
00173 void StreeterFibreGenerator<SPACE_DIM>::PreWriteCalculations(OutputFileHandler& rOutputDirectory)
00174 {
00175 assert(SPACE_DIM == 3);
00176 if (mpGeometryInfo == NULL)
00177 {
00178 EXCEPTION("Files defining the heart surfaces not set");
00179 }
00180
00181
00182 out_stream p_regions_file, p_thickness_file, p_ave_thickness_file;
00183
00184
00185 bool logInfo = PetscTools::AmMaster() & mLogInfo;
00186
00187 if (logInfo)
00188 {
00189 p_regions_file = rOutputDirectory.OpenOutputFile("node_regions.data");
00190 p_thickness_file = rOutputDirectory.OpenOutputFile("wall_thickness.data");
00191 p_ave_thickness_file = rOutputDirectory.OpenOutputFile("averaged_thickness.data");
00192 }
00193
00194
00195 if (fabs(norm_2(mApexToBase)) < DBL_EPSILON)
00196 {
00197 EXCEPTION("Apex to base vector has not been set");
00198 }
00199
00200
00201 unsigned num_nodes = this->mpMesh->GetNumNodes();
00202 for (unsigned node_index=0; node_index<num_nodes; node_index++)
00203 {
00204 double dist_epi, dist_endo;
00205
00206 HeartRegionType node_region = mpGeometryInfo->GetHeartRegion(node_index);
00207
00208 switch(node_region)
00209 {
00210 case HeartGeometryInformation<SPACE_DIM>::LEFT_VENTRICLE_SURFACE:
00211 case HeartGeometryInformation<SPACE_DIM>::LEFT_VENTRICLE_WALL:
00212 dist_epi = mpGeometryInfo->rGetDistanceMapEpicardium()[node_index];
00213 dist_endo = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
00214 break;
00215
00216 case HeartGeometryInformation<SPACE_DIM>::RIGHT_VENTRICLE_SURFACE:
00217 case HeartGeometryInformation<SPACE_DIM>::RIGHT_VENTRICLE_WALL:
00218 dist_epi = mpGeometryInfo->rGetDistanceMapEpicardium()[node_index];
00219 dist_endo = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
00220 break;
00221
00222 case HeartGeometryInformation<SPACE_DIM>::LEFT_SEPTUM:
00223 dist_epi = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
00224 dist_endo = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
00225 break;
00226
00227 case HeartGeometryInformation<SPACE_DIM>::RIGHT_SEPTUM:
00228 dist_epi = mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index];
00229 dist_endo = mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index];
00230 break;
00231
00232 case HeartGeometryInformation<SPACE_DIM>::UNKNOWN:
00233 #define COVERAGE_IGNORE
00234 std::cerr << "Wrong distances node: " << node_index << "\t"
00235 << "Epi " << mpGeometryInfo->rGetDistanceMapEpicardium()[node_index] << "\t"
00236 << "RV " << mpGeometryInfo->rGetDistanceMapRightVentricle()[node_index] << "\t"
00237 << "LV " << mpGeometryInfo->rGetDistanceMapLeftVentricle()[node_index]
00238 << std::endl;
00239
00240
00241 dist_epi = 1;
00242 dist_endo = 0;
00243 break;
00244 #undef COVERAGE_IGNORE
00245
00246 default:
00247 NEVER_REACHED;
00248 }
00249
00250 mWallThickness[node_index] = dist_endo / (dist_endo + dist_epi);
00251
00252 if (std::isnan(mWallThickness[node_index]))
00253 {
00254 #define COVERAGE_IGNORE
00255
00256
00257
00258
00259 mWallThickness[node_index] = 0;
00260 #undef COVERAGE_IGNORE
00261 }
00262
00263 if (logInfo)
00264 {
00265 *p_regions_file << node_region*100 << "\n";
00266 *p_thickness_file << mWallThickness[node_index] << "\n";
00267 }
00268 }
00269
00270
00271
00272
00273 std::vector<double> my_averaged_wall_thickness(num_nodes);
00274 for (unsigned node_index=0; node_index<num_nodes; node_index++)
00275 {
00276 my_averaged_wall_thickness[node_index] = GetAveragedThicknessLocalNode(node_index, mWallThickness);
00277 }
00278
00279
00280 MPI_Allreduce(&my_averaged_wall_thickness[0], &mAveragedWallThickness[0], num_nodes,
00281 MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
00282
00283 if (logInfo)
00284 {
00285 for (unsigned node_index=0; node_index<num_nodes; node_index++)
00286 {
00287 *p_ave_thickness_file << mAveragedWallThickness[node_index] << "\n";
00288 }
00289 }
00290
00291 if (logInfo)
00292 {
00293 p_regions_file->close();
00294 p_thickness_file->close();
00295 p_ave_thickness_file->close();
00296 }
00297 }
00298
00299 template<unsigned SPACE_DIM>
00300 void StreeterFibreGenerator<SPACE_DIM>::Visit(Element<SPACE_DIM, SPACE_DIM>* pElement,
00301 unsigned localElementIndex,
00302 c_vector<double, SPACE_DIM*SPACE_DIM>& rData)
00303 {
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318 c_vector<double, SPACE_DIM> grad_ave_wall_thickness;
00319 c_vector<double, SPACE_DIM+1> elem_nodes_ave_thickness;
00320 double element_averaged_thickness = 0.0;
00321 c_vector<HeartRegionType, SPACE_DIM+1> elem_nodes_region;
00322
00323 for (unsigned local_node_index=0; local_node_index<SPACE_DIM+1; local_node_index++)
00324 {
00325
00326 unsigned global_node_index = pElement->GetNode(local_node_index)->GetIndex();
00327
00328 elem_nodes_ave_thickness[local_node_index] = mAveragedWallThickness[global_node_index];
00329 elem_nodes_region[local_node_index] = mpGeometryInfo->GetHeartRegion(global_node_index);
00330
00331
00332 element_averaged_thickness += mWallThickness[global_node_index];
00333 }
00334
00335 element_averaged_thickness /= SPACE_DIM+1;
00336
00337 c_matrix<double, SPACE_DIM+1, SPACE_DIM> basis_functions( zero_matrix<double>(4u,3u) );
00338 basis_functions(0,0) = basis_functions(0,1) = basis_functions(0,2) = -1.0;
00339 basis_functions(1,0) = basis_functions(2,1) = basis_functions(3,2) = 1.0;
00340
00341 c_matrix<double, SPACE_DIM+1, SPACE_DIM> temp;
00342 c_matrix<double, SPACE_DIM, SPACE_DIM> jacobian, inverse_jacobian;
00343 double jacobian_det;
00344 unsigned element_index = pElement->GetIndex();
00345 this->mpMesh->GetInverseJacobianForElement(element_index, jacobian, jacobian_det, inverse_jacobian);
00346 noalias(temp) = prod (basis_functions, inverse_jacobian);
00347 noalias(grad_ave_wall_thickness) = prod(elem_nodes_ave_thickness, temp);
00348
00349 grad_ave_wall_thickness /= norm_2(grad_ave_wall_thickness);
00350
00351
00352
00353
00354
00355
00356
00357
00358 c_vector<double, SPACE_DIM> fibre_direction = VectorProduct(grad_ave_wall_thickness, mApexToBase);
00359 fibre_direction /= norm_2(fibre_direction);
00360
00361
00362
00363
00364 c_vector<double, SPACE_DIM> longitude_direction = VectorProduct(grad_ave_wall_thickness, fibre_direction);
00365
00366
00367
00368
00369
00370
00371
00372
00373 double alpha = GetFibreMaxAngle(elem_nodes_region) * SmallPow( (1 - 2*element_averaged_thickness), 3 );
00374
00375
00376
00377
00378
00379
00380
00381
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393
00394
00395
00396
00397
00398
00399 c_vector<double, SPACE_DIM> rotated_fibre_direction = fibre_direction*cos(alpha) + longitude_direction*sin(alpha);
00400 c_vector<double, SPACE_DIM> rotated_longitude_direction = -fibre_direction*sin(alpha) + longitude_direction*cos(alpha);
00401
00402
00403
00404
00405
00406 assert( fabs(norm_2(rotated_fibre_direction) - 1) < 100*DBL_EPSILON );
00407 assert( fabs(norm_2(grad_ave_wall_thickness) - 1) < 100*DBL_EPSILON );
00408 assert( fabs(norm_2(rotated_longitude_direction) - 1) < 100*DBL_EPSILON );
00409
00410 assert( fabs(inner_prod(rotated_fibre_direction, grad_ave_wall_thickness)) < 100*DBL_EPSILON );
00411 assert( fabs(inner_prod(rotated_fibre_direction, rotated_longitude_direction)) < 100*DBL_EPSILON);
00412 assert( fabs(inner_prod(grad_ave_wall_thickness, rotated_longitude_direction)) < 100*DBL_EPSILON);
00413
00414
00415
00416
00417
00418
00419
00420
00421 rData[0] = rotated_fibre_direction[0];
00422 rData[1] = rotated_fibre_direction[1];
00423 rData[2] = rotated_fibre_direction[2];
00424 rData[3] = grad_ave_wall_thickness[0];
00425 rData[4] = grad_ave_wall_thickness[1];
00426 rData[5] = grad_ave_wall_thickness[2];
00427 rData[6] = rotated_longitude_direction[0];
00428 rData[7] = rotated_longitude_direction[1];
00429 rData[8] = rotated_longitude_direction[2];
00430 }
00431
00432
00433 template<unsigned SPACE_DIM>
00434 void StreeterFibreGenerator<SPACE_DIM>::SetApexToBase(const c_vector<double, SPACE_DIM>& apexToBase)
00435 {
00436 double norm = norm_2(apexToBase);
00437 if (norm < DBL_EPSILON)
00438 {
00439 EXCEPTION("Apex to base vector should be non-zero");
00440 }
00441 mApexToBase = apexToBase / norm;
00442 }
00443
00444 template<unsigned SPACE_DIM>
00445 void StreeterFibreGenerator<SPACE_DIM>::SetApexToBase(unsigned axis)
00446 {
00447 if (axis >= SPACE_DIM)
00448 {
00449 EXCEPTION("Apex to base coordinate axis was out of range");
00450 }
00451 mApexToBase = zero_vector<double>(SPACE_DIM);
00452 mApexToBase[axis] = 1.0;
00453 }
00454
00455 template<unsigned SPACE_DIM>
00456 void StreeterFibreGenerator<SPACE_DIM>::SetLogInfo(bool logInfo)
00457 {
00458 mLogInfo = logInfo;
00459 }
00460
00462
00464 #define COVERAGE_IGNORE
00465 template class StreeterFibreGenerator<3>;
00466 #undef COVERAGE_IGNORE