Chaste  Build::
VesselNetworkGenerator.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 <sstream>
37 #include <boost/lexical_cast.hpp>
38 #include "RandomNumberGenerator.hpp"
39 #include "SmartPointers.hpp"
40 #include "VoronoiGenerator.hpp"
41 #include "Polygon.hpp"
42 #include "Exception.hpp"
43 #include "VesselNetworkGenerator.hpp"
44 
45 template<unsigned DIM>
47  mReferenceLength(1.e-6 * unit::metres)
48 {
49 }
50 
51 template<unsigned DIM>
53 {
54 }
55 
56 template <unsigned DIM>
57 void VesselNetworkGenerator<DIM>::SetReferenceLengthScale(units::quantity<unit::length> rReferenceLength)
58 {
59  mReferenceLength = rReferenceLength;
60 }
61 
62 template<unsigned DIM>
63 boost::shared_ptr<VesselNetwork<DIM> > VesselNetworkGenerator<DIM>::GenerateOvalNetwork(units::quantity<unit::length> scale_factor,
64  unsigned num_increments,
65  double a_param,
66  double b_param)
67 {
68  // It is 'melon' shaped with one input and output vessel
69  double increment_size = (2.0 * M_PI) / double(num_increments);
70 
71  std::vector<boost::shared_ptr<VesselNode<DIM> > > v1_nodes;
72  v1_nodes.push_back(VesselNode<DIM>::Create(0.0, 0.0, 0.0, mReferenceLength));
73  v1_nodes.push_back(VesselNode<DIM>::Create(increment_size * scale_factor/mReferenceLength, 0.0, 0.0, mReferenceLength));
74  boost::shared_ptr<Vessel<DIM> > p_vessel_1 = Vessel<DIM>::Create(v1_nodes);
75  p_vessel_1->SetId(1);
76 
77  double x_offset = increment_size + std::sqrt(b_param*b_param + a_param*a_param);
78  c_vector<double, 2*DIM> bounds = zero_vector<double>(2*DIM);
79  bounds[0] = 0.0;
80  bounds[4] = 0.0;
81  bounds[5] = 2.0;
82 
83  std::vector<boost::shared_ptr<VesselNode<DIM> > > v2_nodes;
84  double y_max = 0.0;
85  for(unsigned idx=0; idx<= num_increments/2; idx++)
86  {
87  double t = double(idx) * increment_size;
88  double c_value = a_param*a_param*std::cos(2.0 * t) +
89  std::sqrt(pow(b_param,4) - a_param*a_param*pow(std::sin(2.0 * t),2));
90  double x = std::cos(t) * std::sqrt(c_value) + x_offset;
91  double y = std::sin(t) * std::sqrt(c_value);
92  if(y>y_max)
93  {
94  y_max=y;
95  }
96  v2_nodes.push_back(VesselNode<DIM>::Create(x * scale_factor/mReferenceLength, y * scale_factor/mReferenceLength, 0.0, mReferenceLength));
97  }
98 
99  std::vector<boost::shared_ptr<VesselNode<DIM> > > v3_nodes;
100  for(unsigned idx=num_increments/2; idx<= num_increments; idx++)
101  {
102  double t = double(idx) * increment_size;
103  double c_value = a_param*a_param*std::cos(2.0 * t) +
104  std::sqrt(pow(b_param,4) - a_param*a_param*pow(std::sin(2.0 * t),2));
105  double x = std::cos(t) * std::sqrt(c_value) + x_offset;
106  double y = std::sin(t) * std::sqrt(c_value);
107  if(y>y_max)
108  {
109  y_max=y;
110  }
111  v3_nodes.push_back(VesselNode<DIM>::Create(x * scale_factor/mReferenceLength, y * scale_factor/mReferenceLength, 0.0, mReferenceLength));
112  }
113 
114  boost::shared_ptr<Vessel<DIM> > p_vessel_2 = Vessel<DIM>::Create(v2_nodes);
115  boost::shared_ptr<Vessel<DIM> > p_vessel_3 = Vessel<DIM>::Create(v3_nodes);
116  boost::shared_ptr<VesselNetwork<DIM> > p_network = VesselNetwork<DIM>::Create();
117  p_vessel_2->SetId(2);
118  p_vessel_3->SetId(3);
119 
120  std::vector<boost::shared_ptr<VesselNode<DIM> > > v4_nodes;
121  v4_nodes.push_back(VesselNode<DIM>::Create((std::sqrt(a_param*a_param + b_param*b_param) + x_offset) * scale_factor/mReferenceLength, 0.0, 0.0, mReferenceLength));
122  v4_nodes.push_back(VesselNode<DIM>::Create((std::sqrt(a_param*a_param + b_param*b_param) + increment_size + x_offset) * scale_factor/mReferenceLength, 0.0, 0.0, mReferenceLength));
123  bounds[1] = (std::sqrt(a_param*a_param + b_param*b_param) + increment_size + x_offset) * scale_factor/mReferenceLength;
124  bounds[2] = -y_max * scale_factor/mReferenceLength;
125  bounds[3] = y_max * scale_factor/mReferenceLength;
126 
127  boost::shared_ptr<Vessel<DIM> > p_vessel_4 = Vessel<DIM>::Create(v4_nodes);
128  p_vessel_4->SetId(4);
129 
130  p_network->AddVessel(p_vessel_1);
131  p_network->AddVessel(p_vessel_2);
132  p_network->AddVessel(p_vessel_3);
133  p_network->AddVessel(p_vessel_4);
134  p_network->MergeCoincidentNodes(1.e-6);
135 
136  p_network->SetSegmentProperties(p_vessel_1->GetSegments()[0]);
137  p_vessel_1->GetStartNode()->GetFlowProperties()->SetIsInputNode(true);
138  p_vessel_4->GetEndNode()->GetFlowProperties()->SetIsOutputNode(true);
139  return p_network;
140 }
141 
142 template<unsigned DIM>
143 boost::shared_ptr<VesselNetwork<DIM> > VesselNetworkGenerator<DIM>::GenerateParrallelNetwork(boost::shared_ptr<Part<DIM> > domain,
144  units::quantity<unit::per_area> targetDensity,
145  VesselDistribution::Value distributionType,
146  units::quantity<unit::length> exclusionDistance,
147  bool useBbox,
148  std::vector<boost::shared_ptr<DimensionalChastePoint<DIM> > > seeds)
149 {
150  // Get the bounding box of the domain and the volume of the bbox
151  std::vector<units::quantity<unit::length> > bbox = domain->GetBoundingBox();
152  units::quantity<unit::length> delta_x = bbox[1] - bbox[0];
153  units::quantity<unit::length> delta_y = bbox[3] - bbox[2];
154  units::quantity<unit::length> delta_z = 0.0*unit::metres;
155  if(DIM==3)
156  {
157  delta_z = bbox[5] - bbox[4];
158  }
159  if(delta_x == 0.0*unit::metres || delta_y == 0.0*unit::metres)
160  {
161  EXCEPTION("The domain must be at least two-dimensional.");
162  }
163 
164  unsigned num_x = unsigned(units::sqrt(targetDensity) * delta_x);
165  unsigned num_y = unsigned(units::sqrt(targetDensity) * delta_y);
166  if(num_x == 0 || num_y == 0)
167  {
168  EXCEPTION("The domain is not large enough to contain any vessels at the requested density.");
169  }
170  units::quantity<unit::length> spacing_x = delta_x / double(num_x);
171  units::quantity<unit::length> spacing_y = delta_y / double(num_y);
172 
173  // Generate the start and end nodes
174  std::vector<boost::shared_ptr<VesselNode<DIM> > > start_nodes;
175  std::vector<boost::shared_ptr<VesselNode<DIM> > > end_nodes;
176  if(distributionType == VesselDistribution::REGULAR)
177  {
178  for(unsigned idx=0; idx<num_y; idx++)
179  {
180  for(unsigned jdx=0; jdx<num_x; jdx++)
181  {
182  units::quantity<unit::length> location_x = bbox[0] + spacing_x / 2.0 + double(jdx) * spacing_x;
183  units::quantity<unit::length> location_y = bbox[2] + spacing_y / 2.0 + double(idx) * spacing_y;
184  if(DIM==2)
185  {
186  start_nodes.push_back(VesselNode<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, 0.0, mReferenceLength));
187  end_nodes.push_back(VesselNode<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, 0.0, mReferenceLength));
188  }
189  else
190  {
191  start_nodes.push_back(VesselNode<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, bbox[4]/mReferenceLength, mReferenceLength));
192  end_nodes.push_back(VesselNode<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, (bbox[4] + delta_z)/mReferenceLength, mReferenceLength));
193  }
194  }
195  }
196  }
197  else if(distributionType == VesselDistribution::UNIFORM)
198  {
199  unsigned attempts = 0;
200  for(unsigned idx=0; idx<1.e9; idx++)
201  {
202  // Generate with uniform random positions
203  units::quantity<unit::length> location_x = bbox[0] + RandomNumberGenerator::Instance()->ranf() * delta_x;
204  units::quantity<unit::length> location_y = bbox[2] + RandomNumberGenerator::Instance()->ranf() * delta_y;
205 
206  // Get the distance to existing vessels and the boundaries
207  bool free_space = true;
208  bool outside_x = location_x - bbox[0] < exclusionDistance || bbox[1] - location_x < exclusionDistance;
209  bool outside_y = location_y - bbox[2] < exclusionDistance || bbox[3] - location_y < exclusionDistance;
210 
211  if(outside_x || outside_y)
212  {
213  free_space = false;
214  attempts++;
215  }
216  else
217  {
218  for(unsigned kdx=0; kdx<start_nodes.size();kdx++)
219  {
220  double sq_distance = pow((start_nodes[kdx]->rGetLocation().GetLocation(mReferenceLength)[1]- location_y/mReferenceLength),2) +
221  pow((start_nodes[kdx]->rGetLocation().GetLocation(mReferenceLength)[0]- location_x/mReferenceLength),2);
222  if(sq_distance*mReferenceLength*mReferenceLength < (exclusionDistance * exclusionDistance))
223  {
224  free_space = false;
225  attempts++;
226  break;
227  }
228  }
229  }
230 
231  if(free_space)
232  {
233  if(DIM==2)
234  {
235  start_nodes.push_back(VesselNode<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength ,0.0, mReferenceLength));
236  end_nodes.push_back(VesselNode<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, 0.0, mReferenceLength));
237  }
238  else
239  {
240  start_nodes.push_back(VesselNode<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, bbox[4]/mReferenceLength, mReferenceLength));
241  end_nodes.push_back(VesselNode<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, (bbox[4] + delta_z)/mReferenceLength, mReferenceLength));
242  }
243  attempts = 0;
244  }
245  if(start_nodes.size() == num_x * num_y)
246  {
247  break;
248  }
249  if(attempts == 1000)
250  {
251  EXCEPTION("Too many attempts to locate a vessel");
252  }
253  }
254  }
255  else if(distributionType == VesselDistribution::TWO_LAYER)
256  {
257  unsigned attempts = 0;
258 
259  // Uniformly distribute kernels
260  unsigned num_kernels = 8;
261  std::vector<std::vector<units::quantity<unit::length> > > kernel_locations;
262  for(unsigned kdx=0; kdx<num_kernels; kdx++)
263  {
264  std::vector<units::quantity<unit::length> > location;
265  location.push_back(bbox[0] + RandomNumberGenerator::Instance()->ranf() * delta_x);
266  location.push_back(bbox[2] + RandomNumberGenerator::Instance()->ranf() * delta_y);
267  kernel_locations.push_back(location);
268  }
269 
270  // Pick locations randomly from the kernels
271  for(unsigned jdx=0;jdx<1.e9;jdx++)
272  {
273  units::quantity<unit::length> deviation = 100.0* 1.e-6*unit::metres;
274  units::quantity<unit::length> location_x = RandomNumberGenerator::Instance()->NormalRandomDeviate(0.0, deviation/mReferenceLength)*mReferenceLength;
275  units::quantity<unit::length> location_y = RandomNumberGenerator::Instance()->NormalRandomDeviate(0.0, deviation/mReferenceLength)*mReferenceLength;
276  unsigned kernel_index = RandomNumberGenerator::Instance()->randMod(num_kernels);
277  location_x += kernel_locations[kernel_index][0];
278  location_y += kernel_locations[kernel_index][1];
279 
280  // Get the distance to existing vessels and the boundaries
281  bool free_space = true;
282  bool outside_x = location_x - bbox[0] < exclusionDistance || bbox[1] - location_x < exclusionDistance;
283  bool outside_y = location_y - bbox[2] < exclusionDistance || bbox[3] - location_y < exclusionDistance;
284 
285  if(outside_x || outside_y)
286  {
287  free_space = false;
288  attempts++;
289  }
290  else
291  {
292  for(unsigned kdx=0; kdx<start_nodes.size();kdx++)
293  {
294  double sq_distance = pow((start_nodes[kdx]->rGetLocation().GetLocation(mReferenceLength)[1]- location_y/mReferenceLength),2) +
295  pow((start_nodes[kdx]->rGetLocation().GetLocation(mReferenceLength)[0]- location_x/mReferenceLength),2);
296  if(sq_distance*mReferenceLength*mReferenceLength < (exclusionDistance * exclusionDistance))
297  {
298  free_space = false;
299  attempts++;
300  break;
301  }
302  }
303  }
304 
305  if(free_space)
306  {
307  if(DIM==2)
308  {
309  start_nodes.push_back(VesselNode<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, 0.0, mReferenceLength));
310  end_nodes.push_back(VesselNode<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, 0.0, mReferenceLength));
311  }
312  else
313  {
314  start_nodes.push_back(VesselNode<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, bbox[4]/mReferenceLength, mReferenceLength));
315  end_nodes.push_back(VesselNode<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, (bbox[4] + delta_z)/mReferenceLength, mReferenceLength));
316  }
317  attempts = 0;
318  }
319  if(start_nodes.size() == num_x * num_y)
320  {
321  break;
322  }
323  if(attempts == 1000)
324  {
325  EXCEPTION("Too many attempts to locate a vessel");
326  }
327  }
328  }
329 
330  // Set up the vessels
331  std::vector<boost::shared_ptr<Vessel<DIM> > > vessels;
332  for(unsigned idx=0; idx<start_nodes.size(); idx++)
333  {
334  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(start_nodes[idx], end_nodes[idx])));
335  }
336 
337  // Generate and write the network
338  boost::shared_ptr<VesselNetwork<DIM> > p_network = VesselNetwork<DIM>::Create();
339  p_network->AddVessels(vessels);
340 
341  return p_network;
342 }
343 
344 template<unsigned DIM>
345 boost::shared_ptr<VesselNetwork<DIM> > VesselNetworkGenerator<DIM>::Generate3dNetwork(boost::shared_ptr<Part<DIM> > domain,
346  std::vector<units::quantity<unit::per_area> > targetDensity,
347  VesselDistribution::Value distributionType,
348  units::quantity<unit::length> exclusionDistance,
349  bool useBbox,
350  std::vector<boost::shared_ptr<DimensionalChastePoint<DIM> > > seeds)
351 {
352  if(DIM!=3)
353  {
354  EXCEPTION("This generator works in 3-d only.");
355  }
356 
357  if(targetDensity.size()!=3)
358  {
359  EXCEPTION("Target densities must be supplied for each x, y, z dimension.");
360  }
361 
362  // Generate the network
363  boost::shared_ptr<VesselNetwork<DIM> > p_network = VesselNetwork<DIM>::Create();
364 
365  // Get the bounding box of the domain and the volume of the bbox
366  std::vector<units::quantity<unit::length> > bbox = domain->GetBoundingBox();
367  units::quantity<unit::length> delta_x = bbox[1] - bbox[0];
368  units::quantity<unit::length> delta_y = bbox[3] - bbox[2];
369  units::quantity<unit::length> delta_z = bbox[5] - bbox[4];
370 
371  unsigned num_x = unsigned(units::sqrt(targetDensity[0]) * delta_x);
372  unsigned num_y = unsigned(units::sqrt(targetDensity[1]) * delta_y);
373  unsigned num_z = unsigned(units::sqrt(targetDensity[2]) * delta_z);
374  if(num_x == 0 || num_y == 0 || num_z == 0)
375  {
376  EXCEPTION("The domain is not large enough to contain any vessels at the requested density.");
377  }
378  units::quantity<unit::length> spacing_x = delta_x / double(num_x);
379  units::quantity<unit::length> spacing_y = delta_y / double(num_y);
380  units::quantity<unit::length> spacing_z = delta_z / double(num_z);
381 
382  // Generate the nodes
383  std::vector<boost::shared_ptr<VesselNode<DIM> > > nodes;
384  if(distributionType == VesselDistribution::REGULAR)
385  {
386  for(unsigned kdx=0; kdx<num_z+2; kdx++)
387  {
388  for(unsigned idx=0; idx<num_y+2; idx++)
389  {
390  for(unsigned jdx=0; jdx<num_x+2; jdx++)
391  {
392 
393  units::quantity<unit::length> location_x = bbox[0] + spacing_x / 2.0 + double(jdx) * spacing_x -spacing_x;
394  if(jdx==0 || jdx == num_x+1)
395  {
396  location_x -= spacing_x / 2.0;
397  }
398  if(jdx==0)
399  {
400  location_x += spacing_x;
401  }
402 
403  units::quantity<unit::length> location_y = bbox[2] + spacing_y / 2.0 + double(idx) * spacing_y - spacing_y;
404  if(idx==0 || idx == num_y+1)
405  {
406  location_y -= spacing_y / 2.0;
407  }
408  if(idx==0)
409  {
410  location_y += spacing_y;
411  }
412 
413  units::quantity<unit::length> location_z = bbox[4] + spacing_z / 2.0 + double(kdx) * spacing_z - spacing_z;
414  if(kdx==0 || kdx == num_z+1)
415  {
416  location_z -= spacing_z / 2.0;
417  }
418  if(kdx==0)
419  {
420  location_z += spacing_z;
421  }
422 
423  nodes.push_back(VesselNode<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, location_z/mReferenceLength, mReferenceLength));
424  }
425  }
426  }
427 
428  num_x += 2;
429  num_y += 2;
430  num_z += 2;
431 
432  // Set up the vessels
433  std::vector<boost::shared_ptr<Vessel<DIM> > > vessels;
434  for(unsigned kdx=1; kdx<num_z-1; kdx++)
435  {
436  for(unsigned jdx=1; jdx<num_y-1; jdx++)
437  {
438  for(unsigned idx=1; idx<num_x-1; idx++)
439  {
440  unsigned index = idx + num_x*jdx + num_y*num_x*kdx;
441  unsigned index_left = (idx-1) + num_x*jdx + num_y*num_x*kdx;
442  unsigned index_below = idx + num_x*(jdx-1) + num_y*num_x*kdx;
443  unsigned index_front = idx + num_x*(jdx) + num_y*num_x*(kdx-1);
444  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[index_left], nodes[index])));
445  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[index_below], nodes[index])));
446  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[index_front], nodes[index])));
447  if(idx==num_x-2)
448  {
449  unsigned index_right = (idx+1) + num_x*jdx + num_y*num_x*kdx;
450  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[index], nodes[index_right])));
451  }
452  if(jdx==num_y-2)
453  {
454  unsigned index_up = idx + num_x*(jdx+1) + num_y*num_x*kdx;
455  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[index], nodes[index_up])));
456  }
457  if(kdx==num_z-2)
458  {
459  unsigned index_back = idx + num_x*jdx + num_y*num_x*(kdx+1);
460  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[index], nodes[index_back])));
461  }
462  }
463  }
464  }
465  p_network->AddVessels(vessels);
466  }
467 
468  else if(distributionType == VesselDistribution::UNIFORM)
469  {
470  std::vector<boost::shared_ptr<DimensionalChastePoint<DIM> > > seeds;
471  unsigned attempts = 0;
472  for(unsigned idx=0; idx<1.e9; idx++)
473  {
474  // Generate with uniform random positions
475  units::quantity<unit::length> location_x = bbox[0] + RandomNumberGenerator::Instance()->ranf() * delta_x;
476  units::quantity<unit::length> location_y = bbox[2] + RandomNumberGenerator::Instance()->ranf() * delta_y;
477  units::quantity<unit::length> location_z = bbox[4] + RandomNumberGenerator::Instance()->ranf() * delta_z;
478 
479  // Get the distance to existing vessels and the boundaries
480  bool free_space = true;
481  bool outside_x = location_x - bbox[0] < exclusionDistance || bbox[1] - location_x < exclusionDistance;
482  bool outside_y = location_y - bbox[2] < exclusionDistance || bbox[3] - location_y < exclusionDistance;
483  bool outside_z = location_z - bbox[4] < exclusionDistance || bbox[5] - location_z < exclusionDistance;
484 
485  if(outside_x || outside_y || outside_z)
486  {
487  free_space = false;
488  attempts++;
489  }
490  else
491  {
492  for(unsigned kdx=0; kdx<seeds.size();kdx++)
493  {
494  double sq_distance = pow((seeds[kdx]->GetLocation(mReferenceLength)[1]- location_y/mReferenceLength),2) +
495  pow((seeds[kdx]->GetLocation(mReferenceLength)[0]- location_x/mReferenceLength),2) +
496  pow((seeds[kdx]->GetLocation(mReferenceLength)[2]- location_z/mReferenceLength),2);
497  if(sq_distance*mReferenceLength*mReferenceLength < (exclusionDistance * exclusionDistance))
498  {
499  free_space = false;
500  attempts++;
501  break;
502  }
503  }
504  }
505 
506  if(free_space)
507  {
508  seeds.push_back(DimensionalChastePoint<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, location_z/mReferenceLength, mReferenceLength));
509  attempts = 0;
510  }
511  if(seeds.size() == num_x * num_y * num_z)
512  {
513  break;
514  }
515  if(attempts == 1000)
516  {
517  EXCEPTION("Too many attempts to locate a vessel");
518  }
519  }
520 
521  // Generate a network from the resulting voronoi tesselation
522  VoronoiGenerator<DIM> generator;
523  boost::shared_ptr<Part<DIM> > p_tesselation = generator.Generate(domain, seeds, num_x * num_y * num_z);
524  p_network = GenerateFromPart(p_tesselation);
525  }
526  else if(distributionType == VesselDistribution::TWO_LAYER)
527  {
528  std::vector<boost::shared_ptr<DimensionalChastePoint<DIM> > > seeds;
529  unsigned attempts = 0;
530 
531  // Uniformly distribute kernels
532  unsigned num_kernels = 8;
533  std::vector<std::vector<units::quantity<unit::length> > > kernel_locations;
534  for(unsigned kdx=0; kdx<num_kernels; kdx++)
535  {
536  std::vector<units::quantity<unit::length> > location;
537  location.push_back(bbox[0] + RandomNumberGenerator::Instance()->ranf() * delta_x);
538  location.push_back(bbox[2] + RandomNumberGenerator::Instance()->ranf() * delta_y);
539  location.push_back(bbox[3] + RandomNumberGenerator::Instance()->ranf() * delta_z);
540  kernel_locations.push_back(location);
541  }
542 
543  // Pick locations randomly from the kernels
544  for(unsigned jdx=0;jdx<1.e9;jdx++)
545  {
546  units::quantity<unit::length> deviation = 100.0 * 1.e-6 * unit::metres;
547  units::quantity<unit::length> location_x = RandomNumberGenerator::Instance()->NormalRandomDeviate(0.0, deviation/mReferenceLength)*mReferenceLength;
548  units::quantity<unit::length> location_y = RandomNumberGenerator::Instance()->NormalRandomDeviate(0.0, deviation/mReferenceLength)*mReferenceLength;
549  units::quantity<unit::length> location_z = RandomNumberGenerator::Instance()->NormalRandomDeviate(0.0, deviation/mReferenceLength)*mReferenceLength;
550  unsigned kernel_index = RandomNumberGenerator::Instance()->randMod(num_kernels);
551  location_x += kernel_locations[kernel_index][0];
552  location_y += kernel_locations[kernel_index][1];
553  location_z += kernel_locations[kernel_index][2];
554 
555  // Get the distance to existing vessels and the boundaries
556  bool free_space = true;
557  bool outside_x = location_x - bbox[0] < exclusionDistance || bbox[1] - location_x < exclusionDistance;
558  bool outside_y = location_y - bbox[2] < exclusionDistance || bbox[3] - location_y < exclusionDistance;
559  bool outside_z = location_z - bbox[4] < exclusionDistance || bbox[5] - location_z < exclusionDistance;
560 
561  if(outside_x || outside_y || outside_z)
562  {
563  free_space = false;
564  attempts++;
565  }
566  else
567  {
568  for(unsigned kdx=0; kdx<seeds.size();kdx++)
569  {
570  double sq_distance = pow((seeds[kdx]->GetLocation(mReferenceLength)[1]- location_y/mReferenceLength),2) +
571  pow((seeds[kdx]->GetLocation(mReferenceLength)[0]- location_x/mReferenceLength),2) +
572  pow((seeds[kdx]->GetLocation(mReferenceLength)[2]- location_z/mReferenceLength),2);
573  if(sq_distance*mReferenceLength*mReferenceLength < (exclusionDistance * exclusionDistance))
574  {
575  free_space = false;
576  attempts++;
577  break;
578  }
579  }
580  }
581 
582  if(free_space)
583  {
584  seeds.push_back(DimensionalChastePoint<DIM>::Create(location_x/mReferenceLength, location_y/mReferenceLength, location_z/mReferenceLength, mReferenceLength));
585  attempts = 0;
586  }
587  if(seeds.size() == num_x * num_y * num_z)
588  {
589  break;
590  }
591  if(attempts == 1000)
592  {
593  EXCEPTION("Too many attempts to locate a vessel");
594  }
595  }
596 
597  // Generate a network from the resulting voronoi tesselation
598  VoronoiGenerator<DIM> generator;
599  boost::shared_ptr<Part<DIM> > p_tesselation = generator.Generate(domain, seeds, num_x * num_y * num_z);
600  p_network = GenerateFromPart(p_tesselation);
601  }
602  return p_network;
603 }
604 
605 template<unsigned DIM>
607  std::vector<unsigned> numberOfUnits)
608 {
609  // Get unit dimensions
610  std::pair<DimensionalChastePoint<DIM>, DimensionalChastePoint<DIM> > extents = input_unit->GetExtents();
611 
612  // For each number of units in each dimension copy the original vessels and move the copy to the desired location
613  double num_x = 0;
614  double num_y = 0;
615  double num_z = 0;
616  if(numberOfUnits.size() >= 1)
617  {
618  num_x = numberOfUnits[0];
619  if(numberOfUnits.size() >= 2)
620  {
621  num_y = numberOfUnits[1];
622  if(numberOfUnits.size() >= 3 && DIM ==3)
623  {
624  num_z = numberOfUnits[2];
625  }
626  }
627  }
628 
629  // Keep track of the current vessels
630  std::vector<boost::shared_ptr<Vessel<DIM> > > original_vessels = input_unit->GetVessels();
631  for(unsigned idx =0; idx < num_x; idx++)
632  {
633  DimensionalChastePoint<DIM>translation_vector(double(idx+1) * (extents.second.GetLocation(mReferenceLength)[0] - extents.first.GetLocation(mReferenceLength)[0]), 0.0, 0.0, mReferenceLength);
634  std::vector<boost::shared_ptr<Vessel<DIM> > > copied_vessels = input_unit->CopyVessels(original_vessels);
635  input_unit->Translate(translation_vector, copied_vessels);
636  }
637  input_unit->MergeCoincidentNodes();
638 
639  std::vector<boost::shared_ptr<Vessel<DIM> > > x_transform_vessels = input_unit->GetVessels();
640  for(unsigned idx =0; idx < num_y; idx++)
641  {
642  DimensionalChastePoint<DIM>translation_vector(0.0, double(idx+1) * (extents.second.GetLocation(mReferenceLength)[1] - extents.first.GetLocation(mReferenceLength)[1]), 0.0, mReferenceLength);
643  std::vector<boost::shared_ptr<Vessel<DIM> > > copied_vessels = input_unit->CopyVessels(x_transform_vessels);
644  input_unit->Translate(translation_vector, copied_vessels);
645  }
646  input_unit->MergeCoincidentNodes();
647  std::vector<boost::shared_ptr<Vessel<DIM> > > y_transform_vessels = input_unit->GetVessels();
648 
649  for(unsigned idx =0; idx < num_z; idx++)
650  {
651  DimensionalChastePoint<DIM>translation_vector(0.0, 0.0, double(idx+1) * (extents.second.GetLocation(mReferenceLength)[2] - extents.first.GetLocation(mReferenceLength)[2]), mReferenceLength);
652  std::vector<boost::shared_ptr<Vessel<DIM> > > copied_vessels = input_unit->CopyVessels(y_transform_vessels);
653  input_unit->Translate(translation_vector, copied_vessels);
654  }
655  input_unit->MergeCoincidentNodes();
656 }
657 
658 template<unsigned DIM>
659 boost::shared_ptr<VesselNetwork<DIM> > VesselNetworkGenerator<DIM>::GenerateHexagonalUnit(units::quantity<unit::length> vesselLength)
660 {
661  // Generate the nodes
662  double dimensionless_vessel_length = vesselLength/mReferenceLength;
663  std::vector<boost::shared_ptr<VesselNode<DIM> > > nodes;
664  nodes.push_back(VesselNode<DIM>::Create(0.0, 0.0, 0.0, mReferenceLength)); //0
665  nodes.push_back(VesselNode<DIM>::Create(dimensionless_vessel_length, dimensionless_vessel_length, 0.0, mReferenceLength)); //1
666  nodes.push_back(VesselNode<DIM>::Create(0.0, 2.0 * dimensionless_vessel_length, 0.0, mReferenceLength)); //2
667  nodes.push_back(VesselNode<DIM>::Create(2.0 * dimensionless_vessel_length, dimensionless_vessel_length, 0.0, mReferenceLength)); //3
668  nodes.push_back(VesselNode<DIM>::Create(3.0 * dimensionless_vessel_length, 0.0, 0.0, mReferenceLength)); //4
669  nodes.push_back(VesselNode<DIM>::Create(4.0 * dimensionless_vessel_length, 0.0, 0.0, mReferenceLength)); //5
670  nodes.push_back(VesselNode<DIM>::Create(3.0 * dimensionless_vessel_length, 2.0 * dimensionless_vessel_length, 0.0, mReferenceLength)); //6
671  if (DIM == 3)
672  {
673  nodes.push_back(VesselNode<DIM>::Create(0.0, 0.0, 1.5*dimensionless_vessel_length, mReferenceLength)); //7
674  nodes.push_back(VesselNode<DIM>::Create(dimensionless_vessel_length, dimensionless_vessel_length, 1.5*dimensionless_vessel_length, mReferenceLength)); //8
675  nodes.push_back(VesselNode<DIM>::Create(2.0 * dimensionless_vessel_length, dimensionless_vessel_length, 1.5*dimensionless_vessel_length, mReferenceLength)); //9
676  nodes.push_back(VesselNode<DIM>::Create(3.0 * dimensionless_vessel_length, 0.0, 1.5*dimensionless_vessel_length, mReferenceLength)); //9
677  }
678 
679  // Generate the segments and vessels
680  std::vector<boost::shared_ptr<Vessel<DIM> > > vessels;
681  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[0], nodes[1])));
682  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[1], nodes[2])));
683  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[1], nodes[3])));
684  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[3], nodes[4])));
685  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[4], nodes[5])));
686  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[3], nodes[6])));
687 
688  if (DIM == 3)
689  {
690  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[0], nodes[7])));
691  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[1], nodes[8])));
692  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[3], nodes[9])));
693  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[4], nodes[10])));
694  }
695 
696  // Generate the network
697  boost::shared_ptr<VesselNetwork<DIM> > pNetwork(new VesselNetwork<DIM>());
698  pNetwork->AddVessels(vessels);
699 
700  return pNetwork;
701 }
702 
703 template<unsigned DIM>
704 boost::shared_ptr<VesselNetwork<DIM> > VesselNetworkGenerator<DIM>::GenerateBifurcationUnit(units::quantity<unit::length> vesselLength,
705  DimensionalChastePoint<DIM> startPosition)
706 {
707  // Generate the nodes
708  double dimensionless_vessel_length = vesselLength/mReferenceLength;
709  std::vector<boost::shared_ptr<VesselNode<DIM> > > nodes;
710  nodes.push_back(VesselNode<DIM>::Create(0.0, dimensionless_vessel_length, 0.0, mReferenceLength)); //0
711  nodes.push_back(VesselNode<DIM>::Create(dimensionless_vessel_length, dimensionless_vessel_length, 0.0, mReferenceLength)); //1
712  nodes.push_back(VesselNode<DIM>::Create(2.0 * dimensionless_vessel_length, 2.0 * dimensionless_vessel_length, 0.0, mReferenceLength)); //2
713  nodes.push_back(VesselNode<DIM>::Create(2.0 * dimensionless_vessel_length, 0.0, 0.0, mReferenceLength)); //3
714  nodes.push_back(VesselNode<DIM>::Create(3.0 * dimensionless_vessel_length, dimensionless_vessel_length, 0.0, mReferenceLength)); //4
715  nodes.push_back(VesselNode<DIM>::Create(4.0 * dimensionless_vessel_length, dimensionless_vessel_length, 0.0, mReferenceLength)); //5
716 
717  // Generate the segments and vessels
718  std::vector<boost::shared_ptr<Vessel<DIM> > > vessels;
719  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[0], nodes[1])));
720  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[1], nodes[2])));
721  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[1], nodes[3])));
722  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[2], nodes[4])));
723  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[3], nodes[4])));
724  vessels.push_back(Vessel<DIM>::Create(VesselSegment<DIM>::Create(nodes[4], nodes[5])));
725 
726  // Generate the network
727  boost::shared_ptr<VesselNetwork<DIM> > p_network(new VesselNetwork<DIM>());
728  p_network->AddVessels(vessels);
729  p_network->Translate(startPosition);
730  return p_network;
731 }
732 
733 template<unsigned DIM>
734 boost::shared_ptr<VesselNetwork<DIM> > VesselNetworkGenerator<DIM>::GenerateSingleVessel(units::quantity<unit::length> vesselLength,
735  DimensionalChastePoint<DIM> startPosition,
736  unsigned divisions, unsigned axis)
737 {
738  double dimensionless_vessel_length = vesselLength/mReferenceLength;
739  std::vector<boost::shared_ptr<VesselNode<DIM> > > nodes;
740  nodes.push_back(VesselNode<DIM>::Create(startPosition)); //0
741  double interval = dimensionless_vessel_length/double(divisions + 1);
742 
743  for(unsigned idx=0;idx<divisions+1;idx++)
744  {
745  c_vector<double, DIM> dimless_position = startPosition.GetLocation(mReferenceLength);
746  if(DIM==2)
747  {
748  if(axis==0)
749  {
750  nodes.push_back(VesselNode<DIM>::Create(dimless_position[0] + interval*double(idx+1), dimless_position[1], 0.0, mReferenceLength)); //1
751  }
752  else
753  {
754  nodes.push_back(VesselNode<DIM>::Create(dimless_position[0], dimless_position[1] + interval*double(idx+1), 0.0, mReferenceLength)); //1
755  }
756 
757  }
758  else
759  {
760  if(axis==0)
761  {
762  nodes.push_back(VesselNode<DIM>::Create(dimless_position[0] + interval*double(idx+1), dimless_position[1], dimless_position[2], mReferenceLength)); //1
763  }
764  else if(axis==1)
765  {
766  nodes.push_back(VesselNode<DIM>::Create(dimless_position[0], dimless_position[1] + interval*double(idx+1), dimless_position[2], mReferenceLength)); //2
767  }
768  else
769  {
770  nodes.push_back(VesselNode<DIM>::Create(dimless_position[0], dimless_position[1], dimless_position[2] + interval*double(idx+1), mReferenceLength)); //3
771  }
772  }
773  }
774 
775  std::vector<boost::shared_ptr<Vessel<DIM> > > vessels;
776  vessels.push_back(Vessel<DIM>::Create(nodes));
777 
778  // Generate the network
779  boost::shared_ptr<VesselNetwork<DIM> > p_network(new VesselNetwork<DIM>());
780  p_network->AddVessels(vessels);
781  return p_network;
782 }
783 
784 template<unsigned DIM>
785 boost::shared_ptr<VesselNetwork<DIM> > VesselNetworkGenerator<DIM>::GenerateFromPart(boost::shared_ptr<Part<DIM> > part)
786 {
787  boost::shared_ptr<VesselNetwork<DIM> > p_network = VesselNetwork<DIM>::Create();
788 
789  // Get the polygons
790  std::vector<boost::shared_ptr<Polygon<DIM> > > polygons = part->GetPolygons();
791  std::vector<boost::shared_ptr<Vessel<DIM> > > vessels;
792  for (unsigned idx = 0; idx < polygons.size(); idx++)
793  {
794  std::vector<boost::shared_ptr<VesselSegment<DIM> > > segments;
795  std::vector<boost::shared_ptr<DimensionalChastePoint<DIM> > > vertices = polygons[idx]->GetVertices();
796  for (unsigned jdx = 1; jdx < vertices.size(); jdx++)
797  {
798  segments.push_back(VesselSegment<DIM>::Create(VesselNode<DIM>::Create(*vertices[jdx-1]), VesselNode<DIM>::Create(*vertices[jdx])));
799  }
800  vessels.push_back(Vessel<DIM>::Create(segments));
801  }
802  p_network->AddVessels(vessels);
803  p_network->MergeCoincidentNodes();
804  return p_network;
805 }
806 
807 template<unsigned DIM>
808 boost::shared_ptr<VesselNetwork<DIM> > VesselNetworkGenerator<DIM>::GenerateVoronoiNetwork(units::quantity<unit::length> cubeX,
809  units::quantity<unit::length> cubeY,
810  units::quantity<unit::length> cubeZ,
811  unsigned numPoints)
812 {
813  if(cubeZ == 0.0 * unit::metres || DIM==2)
814  {
815  EXCEPTION("This generator only works in 3D");
816  }
817 
818  boost::shared_ptr<Part<DIM> > p_part = Part<DIM>::Create();
819  p_part->AddCuboid(cubeX, cubeY, cubeZ, DimensionalChastePoint<DIM>(0.0, 0.0, 0.0));
820  VoronoiGenerator<DIM> generator;
821  boost::shared_ptr<Part<DIM> > p_tesselation = generator.Generate(p_part, std::vector<boost::shared_ptr<DimensionalChastePoint<DIM> > >(), numPoints);
822  boost::shared_ptr<VesselNetwork<DIM> > p_network = GenerateFromPart(p_tesselation);
823  return p_network;
824 }
825 
826 template<unsigned DIM>
827 void VesselNetworkGenerator<DIM>::MapToSphere(boost::shared_ptr<VesselNetwork<DIM> > pInputUnit,
828  units::quantity<unit::length> radius, units::quantity<unit::length> thickess,
829  double azimuthExtent, double polarExtent)
830 {
831  std::pair<DimensionalChastePoint<DIM>, DimensionalChastePoint<DIM> > extents = pInputUnit->GetExtents();
832  std::vector<boost::shared_ptr<VesselNode<DIM> > > nodes = pInputUnit->GetNodes();
833  for(unsigned idx =0; idx<nodes.size(); idx++)
834  {
835  c_vector<double, DIM> node_loc = nodes[idx]->rGetLocation().GetLocation(mReferenceLength);
836  c_vector<double, DIM> offset = (extents.second - extents.first).GetLocation(mReferenceLength);
837 
838  double x_frac = node_loc[0] / offset[0];
839  double azimuth_angle = x_frac * azimuthExtent;
840 
841  double y_frac = (node_loc[1] + 3.0) / offset[1];
842  double polar_angle = y_frac * polarExtent;
843 
844  double dimless_radius = radius/mReferenceLength;
845  if(offset[2]>0.0)
846  {
847  double z_frac = node_loc[2] / offset[2];
848  dimless_radius = dimless_radius - (thickess/mReferenceLength) * z_frac;
849  }
850  DimensionalChastePoint<DIM>new_position(dimless_radius * std::cos(azimuth_angle) * std::sin(polar_angle),
851  dimless_radius * std::cos(polar_angle),
852  dimless_radius * std::sin(azimuth_angle) * std::sin(polar_angle),
854  nodes[idx]->SetLocation(new_position);
855  }
856 }
857 
858 template<unsigned DIM>
859 boost::shared_ptr<VesselNetwork<DIM> > VesselNetworkGenerator<DIM>::GenerateHexagonalNetwork(units::quantity<unit::length> width,
860  units::quantity<unit::length> height,
861  units::quantity<unit::length> vessel_length)
862 {
863  // Vessels are laid out on a regular grid in a hexagonal pattern.
864  // The repeating unit looks like this:
865  // \_/
866  // / \_
867  // There is an extra set of lines along the top to 'close' the pattern.
868  double grid_length = vessel_length/mReferenceLength;
869  double unit_height = 2.0 * grid_length;
870  double unit_width = 4.0 * grid_length;
871 
872  // The pattern may not reach the extents of the target area.
873  unsigned units_in_x_direction = floor((width/mReferenceLength)/ unit_width);
874  unsigned units_in_y_direction = floor((height/mReferenceLength) / unit_height);
875 
876  // If the number of units is less than 1, just make a single unit in that direction
877  if (units_in_x_direction < 1)
878  {
879  units_in_x_direction = 1;
880  }
881 
882  if (units_in_y_direction < 1)
883  {
884  units_in_x_direction = 1;
885  }
886 
887  // Create vessels by looping over the units, x direction is inside loop.
888  boost::shared_ptr<VesselNetwork<DIM> > pVesselNetwork = VesselNetwork<DIM>::Create();
889  for(unsigned jdx = 0; jdx<units_in_y_direction; jdx++)
890  {
891  for(unsigned idx=0; idx<units_in_x_direction; idx++)
892  {
893  pVesselNetwork->AddVessel(Vessel<DIM>::Create(VesselNode<DIM>::Create(double(idx)*unit_width + 0.0,
894  double(jdx)*unit_height + 0.0, 0, mReferenceLength),
895  VesselNode<DIM>::Create(double(idx)*unit_width + grid_length,
896  double(jdx)*unit_height + grid_length, 0, mReferenceLength)));
897 
898  pVesselNetwork->AddVessel(Vessel<DIM>::Create(VesselNode<DIM>::Create(double(idx)*unit_width + 0.0,
899  double(jdx)*unit_height + 2.0*grid_length, 0, mReferenceLength),
900  VesselNode<DIM>::Create(double(idx)*unit_width + grid_length,
901  double(jdx)*unit_height + grid_length, 0, mReferenceLength)));
902 
903  pVesselNetwork->AddVessel(Vessel<DIM>::Create(VesselNode<DIM>::Create(double(idx)*unit_width + grid_length,
904  double(jdx)*unit_height + grid_length, 0, mReferenceLength),
905  VesselNode<DIM>::Create(double(idx)*unit_width + 2.0*grid_length,
906  double(jdx)*unit_height + grid_length, 0, mReferenceLength)));
907 
908  pVesselNetwork->AddVessel(Vessel<DIM>::Create(VesselNode<DIM>::Create(double(idx)*unit_width + 2.0*grid_length,
909  double(jdx)*unit_height + grid_length, 0, mReferenceLength),
910  VesselNode<DIM>::Create(double(idx)*unit_width + 3.0*grid_length,
911  double(jdx)*unit_height + 2.0*grid_length, 0, mReferenceLength)));
912 
913  pVesselNetwork->AddVessel(Vessel<DIM>::Create(VesselNode<DIM>::Create(double(idx)*unit_width + 2.0*grid_length,
914  double(jdx)*unit_height + grid_length, 0, mReferenceLength),
915  VesselNode<DIM>::Create(double(idx)*unit_width + 3.0*grid_length,
916  double(jdx)*unit_height + 0.0, 0, mReferenceLength)));
917 
918  pVesselNetwork->AddVessel(Vessel<DIM>::Create(VesselNode<DIM>::Create(double(idx)*unit_width + 3.0*grid_length,
919  double(jdx)*unit_height + 0.0, 0, mReferenceLength),
920  VesselNode<DIM>::Create(double(idx)*unit_width + 4.0*grid_length,
921  double(jdx)*unit_height + 0.0, 0, mReferenceLength)));
922  }
923  }
924 
925  // Add an extra line of vessels along the top
926  for(unsigned idx=0; idx<units_in_x_direction; idx++)
927  {
928  pVesselNetwork->AddVessel(Vessel<DIM>::Create(VesselNode<DIM>::Create(double(idx)*unit_width + 3.0*grid_length,
929  double(units_in_y_direction)*unit_height, 0, mReferenceLength),
930  VesselNode<DIM>::Create(double(idx)*unit_width + 4.0 *grid_length,
931  double(units_in_y_direction)*unit_height, 0, mReferenceLength)));
932  }
933 
934 
935  pVesselNetwork->MergeCoincidentNodes();
936  return pVesselNetwork;
937 }
938 
939 //Explicit instantiation
940 template class VesselNetworkGenerator<2> ;
941 template class VesselNetworkGenerator<3> ;
boost::shared_ptr< Part< DIM > > Generate(boost::shared_ptr< Part< DIM > > pPart, std::vector< boost::shared_ptr< DimensionalChastePoint< DIM > > > seeds=std::vector< boost::shared_ptr< DimensionalChastePoint< DIM > > >(), unsigned numSeeds=100)
static boost::shared_ptr< VesselNetwork< DIM > > Create()
void SetReferenceLengthScale(units::quantity< unit::length > rReferenceLength)
static boost::shared_ptr< Vessel< DIM > > Create(boost::shared_ptr< VesselSegment< DIM > > pSegment)
Definition: Vessel.cpp:128
c_vector< double, DIM > GetLocation(units::quantity< unit::length > scale)
static boost::shared_ptr< Part< DIM > > Create()
Definition: Part.cpp:63
units::quantity< unit::length > mReferenceLength
Definition: Part.hpp:60