Chaste  Build::
VesselBasedDiscreteSource.cpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4  All rights reserved.
5 
6  University of Oxford means the Chancellor, Masters and Scholars of the
7  University of Oxford, having an administrative office at Wellington
8  Square, Oxford OX1 2JD, UK.
9 
10  This file is part of Chaste.
11 
12  Redistribution and use in source and binary forms, with or without
13  modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24  AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25  IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26  ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27  LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28  CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29  GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30  HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31  LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32  OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34  */
35 
36 #include "VesselBasedDiscreteSource.hpp"
37 #include "AbstractCellPopulation.hpp"
38 #include "VesselNetwork.hpp"
39 #include "GeometryTools.hpp"
40 #include "Element.hpp"
41 
42 template<unsigned DIM>
44  : DiscreteSource<DIM>(),
45  mVesselPermeability(0.0*unit::metre_per_second),
46  mReferenceConcentration(0.0*unit::mole_per_metre_cubed),
47  mReferenceHaematocrit(1.0)
48 {
49 
50 }
51 
52 template<unsigned DIM>
54 {
55 
56 }
57 
58 template<unsigned DIM>
59 boost::shared_ptr<VesselBasedDiscreteSource<DIM> > VesselBasedDiscreteSource<DIM>::Create()
60 {
61  MAKE_PTR(VesselBasedDiscreteSource<DIM>, pSelf);
62  return pSelf;
63 }
64 
65 template<unsigned DIM>
66 std::vector<units::quantity<unit::concentration_flow_rate> > VesselBasedDiscreteSource<DIM>::GetConstantInUMeshValues()
67 {
68  if(!this->mpMesh)
69  {
70  EXCEPTION("A mesh is required for this type of source");
71  }
72  std::vector<units::quantity<unit::concentration_flow_rate> > values(this->mpMesh->GetNumElements(), 0.0*unit::mole_per_metre_cubed_per_second);
73  std::vector<std::vector<boost::shared_ptr<VesselSegment<DIM> > > > element_segment_map = this->mpMesh->GetElementSegmentMap();
74 
75  for(unsigned idx=0; idx<element_segment_map.size(); idx++)
76  {
77  if(element_segment_map[idx].size()>0)
78  {
79  // Get the element nodal locations and element volume
80  Element<DIM, DIM>* p_element = this->mpMesh->GetElement(idx);
81  std::vector<DimensionalChastePoint<DIM> > element_vertices(4);
82  double determinant = 0.0;
83  c_matrix<double, DIM, DIM> jacobian;
84  p_element->CalculateJacobian(jacobian, determinant);
85  units::quantity<unit::volume> element_volume = p_element->GetVolume(determinant) * units::pow<3>(this->mpMesh->GetReferenceLengthScale());
86  if(p_element->GetNumNodes() != 4)
87  {
88  EXCEPTION("Vessel mesh mapping only supported for linear tetrahedral elements.");
89  }
90  for (unsigned jdx = 0; jdx < 4; jdx++)
91  {
92  element_vertices[jdx] = DimensionalChastePoint<DIM>(p_element->GetNodeLocation(jdx), this->mpMesh->GetReferenceLengthScale());
93  }
94 
95  for (unsigned jdx = 0; jdx < element_segment_map[idx].size(); jdx++)
96  {
97 
98  units::quantity<unit::length> length_in_box = LengthOfLineInTetra<DIM>(element_segment_map[idx][jdx]->GetNode(0)->rGetLocation(),
99  element_segment_map[idx][jdx]->GetNode(1)->rGetLocation(), element_vertices);
100 
101  units::quantity<unit::area> surface_area = 2.0*M_PI*element_segment_map[idx][jdx]->GetRadius()*length_in_box;
102 
103  double haematocrit_Ratio = element_segment_map[idx][jdx]->GetFlowProperties()->GetHaematocrit()/mReferenceHaematocrit;
104  values[idx] += mVesselPermeability * (surface_area/element_volume) * mReferenceConcentration * haematocrit_Ratio;
105  }
106  }
107  }
108 
109  return values;
110 }
111 
112 template<unsigned DIM>
113 std::vector<units::quantity<unit::rate> > VesselBasedDiscreteSource<DIM>::GetLinearInUMeshValues()
114 {
115  if(!this->mpMesh)
116  {
117  EXCEPTION("A mesh is required for this type of source");
118  }
119  std::vector<units::quantity<unit::rate> > values(this->mpMesh->GetNumElements(), 0.0*unit::per_second);
120  std::vector<std::vector<boost::shared_ptr<VesselSegment<DIM> > > > element_segment_map = this->mpMesh->GetElementSegmentMap();
121 
122  for(unsigned idx=0; idx<element_segment_map.size(); idx++)
123  {
124  if(element_segment_map[idx].size()>0)
125  {
126  // Get the element nodal locations and element volume
127  Element<DIM, DIM>* p_element = this->mpMesh->GetElement(idx);
128  std::vector<DimensionalChastePoint<DIM> > element_vertices(4);
129  double determinant = 0.0;
130  c_matrix<double, DIM, DIM> jacobian;
131  p_element->CalculateJacobian(jacobian, determinant);
132  units::quantity<unit::volume> element_volume = p_element->GetVolume(determinant) * units::pow<3>(this->mpMesh->GetReferenceLengthScale());
133  if(p_element->GetNumNodes() != 4)
134  {
135  EXCEPTION("Vessel mesh mapping only supported for linear tetrahedral elements.");
136  }
137  for (unsigned jdx = 0; jdx < 4; jdx++)
138  {
139  element_vertices[jdx] = DimensionalChastePoint<DIM>(p_element->GetNodeLocation(jdx), this->mpMesh->GetReferenceLengthScale());
140  }
141 
142  for (unsigned jdx = 0; jdx < element_segment_map[idx].size(); jdx++)
143  {
144  units::quantity<unit::length> length_in_box = LengthOfLineInTetra<DIM>(element_segment_map[idx][jdx]->GetNode(0)->rGetLocation(),
145  element_segment_map[idx][jdx]->GetNode(1)->rGetLocation(), element_vertices);
146 
147  units::quantity<unit::area> surface_area = 2.0*M_PI*element_segment_map[idx][jdx]->GetRadius()*length_in_box;
148  double haematocrit = element_segment_map[idx][jdx]->GetFlowProperties()->GetHaematocrit();
149  if(haematocrit>0.0)
150  {
151  values[idx] += mVesselPermeability * (surface_area/element_volume);
152  }
153  }
154  }
155  }
156 
157  return values;
158 }
159 
160 template<unsigned DIM>
161 std::vector<units::quantity<unit::concentration_flow_rate> > VesselBasedDiscreteSource<DIM>::GetConstantInURegularGridValues()
162 {
163  std::vector<units::quantity<unit::concentration_flow_rate> > values(this->mpRegularGrid->GetNumberOfPoints(), 0.0*unit::mole_per_metre_cubed_per_second);
164  std::vector<std::vector<boost::shared_ptr<VesselSegment<DIM> > > > point_segment_map = this->mpRegularGrid->GetPointSegmentMap();
165  units::quantity<unit::length> grid_spacing = this->mpRegularGrid->GetSpacing();
166  units::quantity<unit::volume> grid_volume = units::pow<3>(grid_spacing);
167  for(unsigned idx=0; idx<point_segment_map.size(); idx++)
168  {
169  for (unsigned jdx = 0; jdx < point_segment_map[idx].size(); jdx++)
170  {
171  units::quantity<unit::length> length_in_box = LengthOfLineInBox<DIM>(point_segment_map[idx][jdx]->GetNode(0)->rGetLocation(),
172  point_segment_map[idx][jdx]->GetNode(1)->rGetLocation(),
173  this->mpRegularGrid->GetLocationOf1dIndex(idx), grid_spacing);
174 
175  units::quantity<unit::area> surface_area = 2.0*M_PI*point_segment_map[idx][jdx]->GetRadius()*length_in_box;
176 
177  double haematocrit_ratio = point_segment_map[idx][jdx]->GetFlowProperties()->GetHaematocrit()/mReferenceHaematocrit;
178  values[idx] += mVesselPermeability * (surface_area/grid_volume) * mReferenceConcentration * haematocrit_ratio;
179  }
180  }
181  return values;
182 }
183 
184 template<unsigned DIM>
185 std::vector<units::quantity<unit::rate> > VesselBasedDiscreteSource<DIM>::GetLinearInURegularGridValues()
186 {
187  std::vector<units::quantity<unit::rate> > values(this->mpRegularGrid->GetNumberOfPoints(), 0.0*unit::per_second);
188  std::vector<std::vector<boost::shared_ptr<VesselSegment<DIM> > > > point_segment_map = this->mpRegularGrid->GetPointSegmentMap(false);
189  units::quantity<unit::length> grid_spacing = this->mpRegularGrid->GetSpacing();
190  units::quantity<unit::volume> grid_volume = units::pow<3>(grid_spacing);
191  for(unsigned idx=0; idx<point_segment_map.size(); idx++)
192  {
193  for (unsigned jdx = 0; jdx < point_segment_map[idx].size(); jdx++)
194  {
195  units::quantity<unit::length> length_in_box = LengthOfLineInBox<DIM>(point_segment_map[idx][jdx]->GetNode(0)->rGetLocation(),
196  point_segment_map[idx][jdx]->GetNode(1)->rGetLocation(),
197  this->mpRegularGrid->GetLocationOf1dIndex(idx), grid_spacing);
198 
199  units::quantity<unit::area> surface_area = 2.0*M_PI*point_segment_map[idx][jdx]->GetRadius()*length_in_box;
200  double haematocrit = point_segment_map[idx][jdx]->GetFlowProperties()->GetHaematocrit();
201  if(haematocrit>0.0)
202  {
203  values[idx] -= mVesselPermeability * (surface_area/grid_volume);
204  }
205  }
206  }
207  return values;
208 }
209 
210 template<unsigned DIM>
211 void VesselBasedDiscreteSource<DIM>::SetVesselPermeability(units::quantity<unit::membrane_permeability> value)
212 {
213  mVesselPermeability = value;
214 }
215 
216 template<unsigned DIM>
217 void VesselBasedDiscreteSource<DIM>::SetReferenceConcentration(units::quantity<unit::concentration> value)
218 {
219  mReferenceConcentration = value;
220 }
221 
222 template<unsigned DIM>
223 void VesselBasedDiscreteSource<DIM>::SetReferenceHaematocrit(units::quantity<unit::dimensionless> value)
224 {
225  mReferenceHaematocrit = value;
226 }
227 
228 // Explicit instantiation
229 template class VesselBasedDiscreteSource<2>;
230 template class VesselBasedDiscreteSource<3>;
units::quantity< unit::membrane_permeability > mVesselPermeability
void SetVesselPermeability(units::quantity< unit::membrane_permeability > value)
boost::shared_ptr< RegularGrid< DIM > > mpRegularGrid
std::vector< units::quantity< unit::concentration_flow_rate > > GetConstantInURegularGridValues()
void SetReferenceHaematocrit(units::quantity< unit::dimensionless > value)
units::quantity< unit::dimensionless > mReferenceHaematocrit
static boost::shared_ptr< VesselBasedDiscreteSource< DIM > > Create()
units::quantity< unit::concentration > mReferenceConcentration
std::vector< units::quantity< unit::rate > > GetLinearInUMeshValues()
std::vector< units::quantity< unit::concentration_flow_rate > > GetConstantInUMeshValues()
boost::shared_ptr< DiscreteContinuumMesh< DIM, DIM > > mpMesh
std::vector< units::quantity< unit::rate > > GetLinearInURegularGridValues()
void SetReferenceConcentration(units::quantity< unit::concentration > value)