36 #define _BACKWARD_BACKWARD_WARNING_H 1 //Cut out the vtk deprecated warning for now (gcc4.3) 37 #include <vtkDoubleArray.h> 38 #include <vtkPointData.h> 39 #include <vtkSmartPointer.h> 41 #include "VesselSegment.hpp" 42 #include "ChastePoint.hpp" 43 #include "LinearSystem.hpp" 44 #include "ReplicatableVector.hpp" 45 #include "UblasMatrixInclude.hpp" 46 #include "UnitCollection.hpp" 47 #include "RegularGridWriter.hpp" 48 #include "GeometryWriter.hpp" 49 #include "GreensFunctionSolver.hpp" 50 #include "BaseUnits.hpp" 52 template<
unsigned DIM>
58 mSubSegmentCoordinates(),
62 mSegmentConcentration(),
68 mSubsegmentCutoff(1.0*
unit::microns)
73 template<
unsigned DIM>
79 template<
unsigned DIM>
85 template<
unsigned DIM>
97 units::quantity<unit::concentration_flow_rate> sink_rate = this->
mpPde->ComputeConstantInUSourceTerm();
98 units::quantity<unit::volume> sink_volume = units::pow<3>(this->
mpRegularGrid->GetSpacing());
99 mSinkRates = std::vector<units::quantity<unit::molar_flow_rate> >(number_of_sinks, sink_rate * sink_volume);
100 units::quantity<unit::molar_flow_rate> total_sink_rate = std::accumulate(mSinkRates.begin(), mSinkRates.end(), 0.0*unit::mole_per_second);
104 units::quantity<unit::diffusivity> diffusivity = this->
mpPde->ComputeIsotropicDiffusionTerm();
105 std::vector<units::quantity<unit::concentration> > sink_demand_per_subsegment(number_of_subsegments, 0.0*this->
mReferenceConcentration);
106 for (
unsigned idx = 0; idx < number_of_subsegments; idx++)
108 for (
unsigned jdx = 0; jdx < number_of_sinks; jdx++)
110 sink_demand_per_subsegment[idx] += ((*mGvt)[idx][jdx] / diffusivity) * mSinkRates[jdx];
118 double tolerance = 1.e-10;
120 mSourceRates = std::vector<units::quantity<unit::molar_flow_rate> >(number_of_subsegments, 0.0*unit::mole_per_second);
124 units::quantity<unit::amount> reference_amount(1.0*unit::moles);
125 LinearSystem linear_system(number_of_subsegments + 1, number_of_subsegments + 1);
126 linear_system.SetKspType(
"bcgs");
128 for (
unsigned iteration = 0; iteration < 10; iteration++)
130 linear_system.AssembleIntermediateLinearSystem();
131 for (
unsigned i = 0; i < number_of_subsegments; i++)
133 linear_system.SetRhsVectorElement(i, (
mSegmentConcentration[i] - sink_demand_per_subsegment[i])/reference_concentration);
136 linear_system.SetRhsVectorElement(number_of_subsegments, -total_sink_rate*(reference_time/reference_amount));
139 for (
unsigned iter = 0; iter < number_of_subsegments; iter++)
141 for (
unsigned jter = 0; jter < number_of_subsegments; jter++)
143 linear_system.SetMatrixElement(iter, jter, ((*
mGvv)[iter][jter] / diffusivity)*(reference_amount/(reference_time*reference_concentration)));
145 linear_system.SetMatrixElement(number_of_subsegments, iter, 1.0);
146 linear_system.SetMatrixElement(iter, number_of_subsegments, 1.0);
148 linear_system.SetMatrixElement(number_of_subsegments, number_of_subsegments, 0.0);
151 linear_system.AssembleFinalLinearSystem();
152 ReplicatableVector soln_repl(linear_system.Solve());
155 std::vector<units::quantity<unit::molar_flow_rate> > solution_vector(number_of_subsegments + 1);
156 for (
unsigned row = 0; row < number_of_subsegments + 1; row++)
158 (solution_vector)[row] = soln_repl[row]*(reference_amount/reference_time);
162 bool all_in_tolerance =
true;
163 for (
unsigned i = 0; i < number_of_subsegments; i++)
165 units::quantity<unit::molar_flow_rate> diff = units::abs(
mSourceRates[i] - solution_vector[i]);
166 if (diff > tolerance*unit::mole_per_second)
168 all_in_tolerance =
false;
173 for (
unsigned i = 0; i < number_of_subsegments; i++)
177 g0 = solution_vector[number_of_subsegments]*reference_concentration*(reference_time/reference_amount);
179 if (all_in_tolerance)
187 std::cout <<
"Did not converge\n";
193 for (
unsigned i = 0; i < number_of_sinks; i++)
196 for (
unsigned j = 0; j < number_of_sinks; j++)
198 this->mConcentrations[i] += (*mGtt)[i][j] * mSinkRates[j] / diffusivity;
201 for (
unsigned j = 0; j < number_of_subsegments; j++)
203 this->mConcentrations[i] += (*mGtv)[i][j] *
mSourceRates[j] / diffusivity;
205 this->mConcentrations[i] += g0;
208 std::map<std::string, std::vector<units::quantity<unit::concentration> > > segmentPointData;
218 template<
unsigned DIM>
224 std::vector<boost::shared_ptr<Vessel<DIM> > > vessels = this->
mpNetwork->GetVessels();
225 typename std::vector<boost::shared_ptr<Vessel<DIM> > >::iterator vessel_iter;
226 typename std::vector<boost::shared_ptr<VesselSegment<DIM> > >::iterator segment_iter;
231 for (vessel_iter = vessels.begin(); vessel_iter != vessels.end(); vessel_iter++)
233 std::vector<boost::shared_ptr<VesselSegment<DIM> > > segments = (*vessel_iter)->GetSegments();
234 for (segment_iter = segments.begin(); segment_iter != segments.end(); segment_iter++)
236 units::quantity<unit::length> segment_length = (*segment_iter)->GetLength();
239 if (segment_length < 1.01 * max_subsegment_length)
251 double subsegment_length_factor = segment_length / max_subsegment_length;
252 unsigned num_subsegments = std::floor(subsegment_length_factor) + 1;
253 units::quantity<unit::length> subsegment_length = segment_length / double(num_subsegments);
255 for (
unsigned i = 0; i < num_subsegments; i++)
266 template<
unsigned DIM>
269 unsigned num_points = this->
mpRegularGrid->GetNumberOfPoints();
272 for(
unsigned idx=0; idx<num_points; idx++)
275 mSinkPointMap[idx] = idx;
279 template<
unsigned DIM>
302 template<
unsigned DIM>
305 typedef boost::multi_array<units::quantity<unit::per_length>, 2>::index index;
307 double coefficient = 1.0 / (4.0 * M_PI);
309 boost::shared_ptr<boost::multi_array<units::quantity<unit::per_length>, 2> > p_interaction_matrix(
new boost::multi_array<units::quantity<unit::per_length>, 2>(boost::extents[num_sub_segments][num_sub_segments]));
310 for (index iter = 0; iter < num_sub_segments; iter++)
312 for (index iter2 = 0; iter2 < num_sub_segments; iter2++)
317 units::quantity<unit::per_length> term;
322 double green_correction = 0.6 * std::exp(-0.45 * max_segment_length /radius);
327 term = (1.298 / (1.0 + 0.297 * pow(max_segment_length /radius, 0.838))- green_correction * pow(distance / radius, 2))*coefficient/radius;
331 term = coefficient / distance;
333 (*p_interaction_matrix)[iter][iter2] = term;
334 (*p_interaction_matrix)[iter2][iter] = term;
338 return p_interaction_matrix;
341 template<
unsigned DIM>
344 typedef boost::multi_array<units::quantity<unit::per_length>, 2>::index index;
346 double coefficient = 1.0 / (4.0 * M_PI);
347 units::quantity<unit::volume> tissue_point_volume = units::pow<3>(this->
mpRegularGrid->GetSpacing());
348 units::quantity<unit::length> equivalent_tissue_point_radius = units::root<3>(tissue_point_volume * 0.75 / M_PI);
349 boost::shared_ptr<boost::multi_array<units::quantity<unit::per_length>, 2 > > p_interaction_matrix(
new boost::multi_array<units::quantity<unit::per_length>, 2>(boost::extents[num_points][num_points]));
350 for (index iter = 0; iter < num_points; iter++)
352 for (index iter2 = 0; iter2 < num_points; iter2++)
357 (*p_interaction_matrix)[iter][iter2] = coefficient / distance;
358 (*p_interaction_matrix)[iter2][iter] = coefficient / distance;
360 else if (iter == iter2)
362 (*p_interaction_matrix)[iter][iter2] = 1.2 * coefficient / equivalent_tissue_point_radius;
366 return p_interaction_matrix;
369 template<
unsigned DIM>
372 typedef boost::multi_array<units::quantity<unit::per_length>, 2>::index index;
376 units::quantity<unit::volume> tissue_point_volume = units::pow<3>(this->
mpRegularGrid->GetSpacing());
377 units::quantity<unit::length> equivalent_tissue_point_radius = units::root<3>(tissue_point_volume * 0.75 / M_PI);
378 double coefficient = 1.0 / (4.0 * M_PI);
380 boost::shared_ptr<boost::multi_array<units::quantity<unit::per_length>, 2> > p_interaction_matrix(
new boost::multi_array<units::quantity<unit::per_length>, 2>(boost::extents[num_sinks][num_subsegments]));
381 for (index iter = 0; iter < num_sinks; iter++)
383 for (index iter2 = 0; iter2 < num_subsegments; iter2++)
386 units::quantity<unit::per_length> term;
387 if (distance <= equivalent_tissue_point_radius)
389 term = coefficient * (1.5 - 0.5 * (pow(distance / equivalent_tissue_point_radius, 2))) / equivalent_tissue_point_radius;
393 term = coefficient / distance;
395 (*p_interaction_matrix)[iter][iter2] = term;
398 return p_interaction_matrix;
401 template<
unsigned DIM>
404 typedef boost::multi_array<units::quantity<unit::per_length>, 2>::index index;
407 double coefficient = 1.0 / (4.0 * M_PI);
409 units::quantity<unit::volume> tissue_point_volume = units::pow<3>(this->
mpRegularGrid->GetSpacing());
410 units::quantity<unit::length> equivalent_tissue_point_radius = units::root<3>(tissue_point_volume * 0.75 / M_PI);
412 boost::shared_ptr<boost::multi_array<units::quantity<unit::per_length>, 2> > p_interaction_matrix(
new boost::multi_array<units::quantity<unit::per_length>, 2>(boost::extents[num_subsegments][num_sinks]));
413 for (index iter = 0; iter < num_subsegments; iter++)
415 for (index iter2 = 0; iter2 < num_sinks; iter2++)
418 units::quantity<unit::per_length> term;
419 if (distance <= equivalent_tissue_point_radius)
421 term = coefficient * (1.5 - 0.5 * (pow(distance / equivalent_tissue_point_radius, 2)))/ equivalent_tissue_point_radius;
425 term = coefficient / distance;
427 (*p_interaction_matrix)[iter][iter2] = term;
430 return p_interaction_matrix;
433 template<
unsigned DIM>
443 vtkSmartPointer<vtkPolyData> pPolyData = vtkSmartPointer<vtkPolyData>::New();
444 vtkSmartPointer<vtkPoints> pPoints = vtkSmartPointer<vtkPoints>::New();
448 pPoints->InsertNextPoint(location[0], location[1], location[2]);
450 pPolyData->SetPoints(pPoints);
453 std::map<std::string, std::vector<units::quantity<unit::concentration> > >::iterator segment_iter;
454 for (segment_iter = segmentPointData.begin(); segment_iter != segmentPointData.end(); ++segment_iter)
456 vtkSmartPointer<vtkDoubleArray> pInfo = vtkSmartPointer<vtkDoubleArray>::New();
457 pInfo->SetNumberOfComponents(1);
459 pInfo->SetName(segment_iter->first.c_str());
463 pInfo->SetValue(i, segmentPointData[segment_iter->first][i]/this->mReferenceConcentration);
465 pPolyData->GetPointData()->AddArray(pInfo);
470 geometry_writer.
SetInput(pPolyData);
471 geometry_writer.
Write();
boost::shared_ptr< boost::multi_array< units::quantity< unit::per_length >, 2 > > mGtt
std::vector< units::quantity< unit::length > > mSubSegmentLengths
boost::shared_ptr< AbstractDiscreteContinuumLinearEllipticPde< DIM, DIM > > mpPde
std::vector< units::quantity< unit::molar_flow_rate > > mSourceRates
units::quantity< unit::length > mSubsegmentCutoff
void SetSubSegmentCutoff(units::quantity< unit::length > value)
std::vector< DimensionalChastePoint< DIM > > mSubSegmentCoordinates
boost::shared_ptr< boost::multi_array< units::quantity< unit::per_length >, 2 > > GetTissueVesselInteractionMatrix()
void SetFilename(const std::string &rFilename)
boost::shared_ptr< boost::multi_array< units::quantity< unit::per_length >, 2 > > mGvt
boost::shared_ptr< RegularGrid< DIM > > mpRegularGrid
std::vector< units::quantity< unit::concentration > > mSegmentConcentration
units::quantity< unit::time > GetReferenceTimeScale()
boost::shared_ptr< OutputFileHandler > mpOutputFileHandler
std::vector< units::quantity< unit::molar_flow_rate > > mSinkRates
std::map< unsigned, boost::shared_ptr< VesselSegment< DIM > > > mSegmentPointMap
void SetInput(vtkSmartPointer< vtkPolyData > pSurface)
std::vector< DimensionalChastePoint< DIM > > mSinkCoordinates
std::vector< units::quantity< unit::concentration > > mConcentrations
boost::shared_ptr< boost::multi_array< units::quantity< unit::per_length >, 2 > > mGtv
vtkSmartPointer< vtkImageData > mpVtkSolution
void GenerateSubSegments()
virtual void UpdateSolution(std::vector< double > &rData)
boost::shared_ptr< boost::multi_array< units::quantity< unit::per_length >, 2 > > mGvv
static BaseUnits * Instance()
void UpdateGreensFunctionMatrices(bool updateGtt=0, bool updateGvv=0, bool updateGtv=0, bool updateGvt=0)
boost::shared_ptr< VesselNetwork< DIM > > mpNetwork
void GenerateTissuePoints()
boost::shared_ptr< boost::multi_array< units::quantity< unit::per_length >, 2 > > GetVesselVesselInteractionMatrix()
units::quantity< unit::concentration > mReferenceConcentration
boost::shared_ptr< boost::multi_array< units::quantity< unit::per_length >, 2 > > GetVesselTissueInteractionMatrix()
void SetImage(vtkSmartPointer< vtkImageData > pImage)
boost::shared_ptr< boost::multi_array< units::quantity< unit::per_length >, 2 > > GetTissueTissueInteractionMatrix()
void SetFileName(const std::string &rFileName)
void WriteSolution(std::map< std::string, std::vector< units::quantity< unit::concentration > > > &segmentPointData)
std::vector< unsigned > mSinkPointMap