Chaste  Build::
GeometryToolsImpl.hpp
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 #define _BACKWARD_BACKWARD_WARNING_H 1 //Cut out the vtk deprecated warning for now (gcc4.3)
37 #include <vtkBox.h>
38 #include <vtkTetra.h>
39 #include <vtkPoints.h>
40 #include <vtkSmartPointer.h>
41 #include <vtkUnstructuredGrid.h>
42 #include <vtkCellLocator.h>
43 #include "Exception.hpp"
44 #include "GeometryTools.hpp"
45 #include "UblasVectorInclude.hpp"
46 #include "UblasIncludes.hpp"
47 #include "UblasCustomFunctions.hpp"
48 
49 template<unsigned DIM>
50 units::quantity<unit::length> GetDistance(const DimensionalChastePoint<DIM>& rLocation1,
51  const DimensionalChastePoint<DIM>& rLocation2)
52 {
53  units::quantity<unit::length> reference_length = rLocation1.GetReferenceLengthScale();
54  return norm_2(rLocation2.GetLocation(reference_length) - rLocation1.GetLocation(reference_length))*reference_length;
55 }
56 
57 template<unsigned DIM>
58 units::quantity<unit::area> GetDotProduct(const DimensionalChastePoint<DIM>& rLocation1,
59  const DimensionalChastePoint<DIM>& rLocation2)
60 {
61  units::quantity<unit::length> reference_length = rLocation1.GetReferenceLengthScale();
62  return inner_prod(rLocation2.GetLocation(reference_length), rLocation1.GetLocation(reference_length))*reference_length*reference_length;
63 }
64 
65 template<unsigned DIM>
66 units::quantity<unit::length> GetDotProduct(const DimensionalChastePoint<DIM>& rLocation1,
67  const c_vector<double, DIM>& rLocation2)
68 {
69  units::quantity<unit::length> reference_length_1 = rLocation1.GetReferenceLengthScale();
70  return inner_prod(rLocation2 , rLocation1.GetLocation(reference_length_1))*reference_length_1;
71 }
72 
73 template<unsigned DIM>
74 DimensionalChastePoint<DIM> GetPointProjectionOnLineSegment(const DimensionalChastePoint<DIM>& rStartLocation,
75  const DimensionalChastePoint<DIM>& rEndLocation,
76  const DimensionalChastePoint<DIM>& rProbeLocation,
77  bool projectToEnds,
78  bool checkDimensions)
79 {
80  DimensionalChastePoint<DIM> segment_vector = rEndLocation - rStartLocation;
81  DimensionalChastePoint<DIM> point_vector = rProbeLocation - rStartLocation;
82  units::quantity<unit::area> dp_segment_point = GetDotProduct(segment_vector, point_vector);
83  units::quantity<unit::area> dp_segment_segment = GetDotProduct(segment_vector, segment_vector);
84 
85  if (dp_segment_point <= 0.0*unit::metres*unit::metres || dp_segment_segment <= dp_segment_point)
86  {
87  if(!projectToEnds)
88  {
89  EXCEPTION("Projection of point is outside segment.");
90  }
91  else
92  {
93  units::quantity<unit::length> dist1 = (rStartLocation - rEndLocation).GetNorm2();
94  units::quantity<unit::length> dist2 = (rEndLocation - rProbeLocation).GetNorm2();
95  if(dist1 <= dist2)
96  {
97  return rStartLocation;
98  }
99  else
100  {
101  return rEndLocation;
102  }
103  }
104  }
105  // Point projection is inside segment, get distance to point projection
106  return rStartLocation + segment_vector*(dp_segment_point / dp_segment_segment);
107 }
108 
109 template<unsigned DIM>
110 units::quantity<unit::length> GetDistanceToLineSegment(const DimensionalChastePoint<DIM>& rStartLocation,
111  const DimensionalChastePoint<DIM>& rEndLocation,
112  const DimensionalChastePoint<DIM>& rProbeLocation)
113 {
114  DimensionalChastePoint<DIM> segment_vector = rEndLocation - rStartLocation;
115  units::quantity<unit::area> dp_segment_point = GetDotProduct(segment_vector, rProbeLocation - rStartLocation);
116  // Point projection is outside segment, return node0 distance
117  if (dp_segment_point <= 0.0*unit::metres*unit::metres)
118  {
119  return rStartLocation.GetDistance(rProbeLocation);
120  }
121 
122  units::quantity<unit::area> dp_segment_segment = GetDotProduct(segment_vector, segment_vector);
123  // Point projection is outside segment, return node1 distance
124  if (dp_segment_segment <= dp_segment_point)
125  {
126  return rEndLocation.GetDistance(rProbeLocation);
127  }
128 
129  // Point projection is inside segment, get distance to point projection
130  double projection_ratio = dp_segment_point / dp_segment_segment;
131  DimensionalChastePoint<DIM> projected_location = rStartLocation + segment_vector*projection_ratio - rProbeLocation;
132  return projected_location.GetNorm2();
133 }
134 
135 template<unsigned DIM>
136 std::vector<DimensionalChastePoint<DIM> > GetProbeLocationsExternalPoint(DimensionalChastePoint<DIM> rCentrePoint,
137  units::quantity<unit::length> probeLength)
138 {
139  unsigned num_probes = (2*DIM)+1;
140  units::quantity<unit::length> length_scale = rCentrePoint.GetReferenceLengthScale();
141  std::vector<DimensionalChastePoint<DIM> > probe_locations(num_probes, DimensionalChastePoint<DIM>(0.0, 0.0, 0.0, length_scale));
142 
143  double normalized_probe_length = probeLength/length_scale;
144  probe_locations[0] = rCentrePoint;
145  probe_locations[1] = rCentrePoint + DimensionalChastePoint<DIM>(normalized_probe_length * unit_vector<double>(DIM,0), length_scale);
146  probe_locations[2] = rCentrePoint + DimensionalChastePoint<DIM>(normalized_probe_length * unit_vector<double>(DIM,0), length_scale);
147  probe_locations[3] = rCentrePoint + DimensionalChastePoint<DIM>(normalized_probe_length * unit_vector<double>(DIM,1), length_scale);
148  probe_locations[4] = rCentrePoint + DimensionalChastePoint<DIM>(normalized_probe_length * unit_vector<double>(DIM,1), length_scale);
149  if(DIM==3)
150  {
151  probe_locations[5] = rCentrePoint + DimensionalChastePoint<DIM>(normalized_probe_length * unit_vector<double>(DIM,2), length_scale);
152  probe_locations[6] = rCentrePoint + DimensionalChastePoint<DIM>(normalized_probe_length * unit_vector<double>(DIM,2), length_scale);
153  }
154  return probe_locations;
155 }
156 
157 template<unsigned DIM>
158 std::vector<DimensionalChastePoint<DIM> > GetProbeLocationsInternalPoint(DimensionalChastePoint<DIM> rInitialDirection,
159  DimensionalChastePoint<DIM> rCentralPoint,
160  DimensionalChastePoint<DIM> rRotationAxis,
161  units::quantity<unit::length> probeLength,
162  units::quantity<unit::plane_angle> angle)
163 {
164  unsigned num_probes = 2*DIM - 1;
165  units::quantity<unit::length> length_scale = rCentralPoint.GetReferenceLengthScale();
166  std::vector<DimensionalChastePoint<DIM> > probe_locations(num_probes, DimensionalChastePoint<DIM>(0.0, 0.0, 0.0, length_scale));
167  probe_locations[0] = rCentralPoint;
168  double normalized_probe_length = probeLength/length_scale;
169 
170  if(DIM==3)
171  {
172  c_vector<double, DIM> unit_axis = rRotationAxis.GetUnitVector();
173  c_vector<double, DIM> new_direction = RotateAboutAxis<DIM>(rInitialDirection.GetUnitVector(), unit_axis, angle);
174  new_direction /= norm_2(new_direction);
175  probe_locations[1] = rCentralPoint + DimensionalChastePoint<DIM>(new_direction*normalized_probe_length, length_scale);
176 
177  c_vector<double, DIM> new_direction_r1 = RotateAboutAxis<DIM>(new_direction, unit_axis, M_PI*unit::radians);
178  new_direction_r1 /= norm_2(new_direction_r1);
179  probe_locations[2] = rCentralPoint + DimensionalChastePoint<DIM>(new_direction_r1*normalized_probe_length, length_scale);
180 
181  c_vector<double, DIM> new_direction_r2 = RotateAboutAxis<DIM>(new_direction, unit_axis, M_PI/2.0*unit::radians);
182  new_direction_r2 /= norm_2(new_direction_r2);
183  probe_locations[3] = rCentralPoint + DimensionalChastePoint<DIM>(new_direction_r2*normalized_probe_length, length_scale);
184 
185  c_vector<double, DIM> new_direction_r3 = RotateAboutAxis<DIM>(new_direction, unit_axis, 3.0*M_PI/2.0*unit::radians);
186  new_direction_r3 /= norm_2(new_direction_r3);
187  probe_locations[4] = rCentralPoint + DimensionalChastePoint<DIM>(new_direction_r3*normalized_probe_length, length_scale);
188  }
189  else
190  {
191  DimensionalChastePoint<DIM> new_direction(normalized_probe_length*rInitialDirection.GetUnitVector(), length_scale);
192  probe_locations[1] = rCentralPoint + new_direction;
193  probe_locations[2] = rCentralPoint - new_direction;
194  }
195 
196  return probe_locations;
197 }
198 
199 template<unsigned DIM>
200 bool IsPointInCone(const DimensionalChastePoint<DIM>& rPoint,
201  const DimensionalChastePoint<DIM>& rApex,
202  const DimensionalChastePoint<DIM>& rBase,
203  double aperture)
204 {
205  DimensionalChastePoint<DIM> apex_to_point = rApex - rPoint;
206  DimensionalChastePoint<DIM> apex_to_base = rApex - rBase;
207  units::quantity<unit::length> dist_apex_base = apex_to_base.GetNorm2();
208  units::quantity<unit::area> dp_point_base = GetDotProduct(apex_to_point, apex_to_base);
209  bool in_infinite_cone = dp_point_base / (apex_to_point.GetNorm2() * dist_apex_base) > std::cos(aperture/2.0);
210  if(!in_infinite_cone)
211  {
212  return false;
213  }
214  return dp_point_base / dist_apex_base < dist_apex_base;
215 }
216 
217 template<unsigned DIM>
218 bool IsPointInBox(const DimensionalChastePoint<DIM>& rPoint,
219  const DimensionalChastePoint<DIM>& rLocation, units::quantity<unit::length> spacing)
220 {
221  bool point_in_box = false;
222  units::quantity<unit::length> point_length_scale = rPoint.GetReferenceLengthScale();
223  c_vector<double, DIM> location_in_point_scale = rLocation.GetLocation(point_length_scale);
224  c_vector<double, DIM> dimensionless_point = rPoint.GetLocation(point_length_scale);
225  double dimensionless_spacing = spacing/point_length_scale;
226 
227  bool inside_left = dimensionless_point[0] >= location_in_point_scale[0] -dimensionless_spacing/2.0;
228  bool inside_right = dimensionless_point[0] <= location_in_point_scale[0] + dimensionless_spacing/2.0;
229  if(inside_left && inside_right)
230  {
231  if(dimensionless_point[1] >= location_in_point_scale[1] -dimensionless_spacing/2.0 && dimensionless_point[1] <=
232  location_in_point_scale[1] + dimensionless_spacing/2.0)
233  {
234  if(DIM==3)
235  {
236  if(dimensionless_point[2] >= location_in_point_scale[2] -dimensionless_spacing/2.0 && dimensionless_point[2] <=
237  location_in_point_scale[2] + dimensionless_spacing/2.0)
238  {
239  return true;
240  }
241  }
242  else
243  {
244  return true;
245  }
246  }
247  }
248  return point_in_box;
249 }
250 
251 template<unsigned DIM>
252 bool IsPointInTetra(const DimensionalChastePoint<DIM>& rPoint, const std::vector<DimensionalChastePoint<DIM> >& locations)
253 {
254  units::quantity<unit::length> scale_factor = rPoint.GetReferenceLengthScale();
255  vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints> :: New();
256  c_vector<double, DIM> loc0 = locations[0].GetLocation(scale_factor);
257  c_vector<double, DIM> loc1 = locations[1].GetLocation(scale_factor);
258  c_vector<double, DIM> loc2 = locations[2].GetLocation(scale_factor);
259  c_vector<double, DIM> loc3 = locations[3].GetLocation(scale_factor);
260  if(DIM==3)
261  {
262  points->InsertNextPoint(&loc0[0]);
263  points->InsertNextPoint(&loc1[0]);
264  points->InsertNextPoint(&loc2[0]);
265  points->InsertNextPoint(&loc3[0]);
266  }
267  else
268  {
269  EXCEPTION("PointInTetra only works in 3d");
270  }
271 
272  vtkSmartPointer<vtkUnstructuredGrid> p_grid = vtkSmartPointer<vtkUnstructuredGrid>::New();
273  p_grid->SetPoints(points);
274  vtkIdType ptIds[] = {0, 1, 2, 3};
275  p_grid->InsertNextCell(VTK_TETRA, 4, ptIds);
276 
277  vtkSmartPointer<vtkCellLocator> p_locator = vtkSmartPointer<vtkCellLocator>::New();
278  p_locator->SetDataSet(p_grid);
279  p_locator->Update();
280  c_vector<double, DIM> probe_location = rPoint.GetLocation(scale_factor);
281  int in_tetra = p_locator->FindCell(&probe_location[0]);
282  if(in_tetra == -1)
283  {
284  return false;
285  }
286  else
287  {
288  return true;
289  }
290 }
291 
292 template<unsigned DIM>
293 units::quantity<unit::length> LengthOfLineInBox(const DimensionalChastePoint<DIM>& rStartPoint,
294  const DimensionalChastePoint<DIM>& rEndPoint,
295  const DimensionalChastePoint<DIM>& rLocation, units::quantity<unit::length> spacing)
296 {
297  // If the line is fully in the box return its length
298  bool point1_in_box = IsPointInBox<DIM>(rStartPoint, rLocation, spacing);
299  bool point2_in_box = IsPointInBox<DIM>(rEndPoint, rLocation, spacing);
300  if(point1_in_box && point2_in_box)
301  {
302  return GetDistance(rEndPoint, rStartPoint);
303  }
304  else
305  {
306  units::quantity<unit::length> scale_factor = rLocation.GetReferenceLengthScale();
307  double dimensionless_spacing = spacing/scale_factor;
308 
309  c_vector<double,6> dimensionless_bounds;
310  c_vector<double, DIM> dimensionless_location = rLocation.GetLocation(scale_factor);
311  dimensionless_bounds[0] = dimensionless_location[0] - dimensionless_spacing/2.0;
312  dimensionless_bounds[1] = dimensionless_location[0] + dimensionless_spacing/2.0;
313  dimensionless_bounds[2] = dimensionless_location[1] - dimensionless_spacing/2.0;
314  dimensionless_bounds[3] = dimensionless_location[1] + dimensionless_spacing/2.0;
315  if(DIM==3)
316  {
317  dimensionless_bounds[4] = dimensionless_location[2] - dimensionless_spacing/2.0;
318  dimensionless_bounds[5] = dimensionless_location[2] + dimensionless_spacing/2.0;
319  }
320  else
321  {
322  dimensionless_bounds[4] = 0.0;
323  dimensionless_bounds[5] = 0.0;
324  }
325 
326  double t1;
327  double t2;
328  int plane1;
329  int plane2;
330  c_vector<double,DIM> intercept_1;
331  c_vector<double,DIM> intercept_2;
332 
333  c_vector<double,3> dimensionless_start = rStartPoint.GetLocation(scale_factor);
334  c_vector<double,3> dimensionless_end = rEndPoint.GetLocation(scale_factor);
335  int in_box = vtkBox::IntersectWithLine(&dimensionless_bounds[0], &dimensionless_start[0], &dimensionless_end[0],
336  t1, t2, &intercept_1[0], &intercept_2[0], plane1, plane2);
337 
338  if(point1_in_box)
339  {
340  return norm_2(intercept_2 - dimensionless_start)*scale_factor;
341  }
342 
343  if(point2_in_box)
344  {
345  return norm_2(intercept_1 - dimensionless_end)*scale_factor;
346  }
347 
348  if(in_box)
349  {
350  return norm_2(intercept_2 - intercept_1)*scale_factor;
351  }
352  else
353  {
354  return 0.0*scale_factor;
355  }
356  }
357 }
358 
359 template<unsigned DIM>
360 units::quantity<unit::length> LengthOfLineInTetra(const DimensionalChastePoint<DIM>& rStartPoint,
361  const DimensionalChastePoint<DIM>& rEndPoint,
362  const std::vector<DimensionalChastePoint<DIM> >& locations)
363 {
364  bool point1_in_tetra = IsPointInTetra<DIM>(rStartPoint, locations);
365  bool point2_in_tetra = IsPointInTetra<DIM>(rEndPoint, locations);
366 
367  if(point1_in_tetra && point2_in_tetra)
368  {
369  return GetDistance(rEndPoint, rStartPoint);
370  }
371  else
372  {
373  int line_crosses;
374 
375  units::quantity<unit::length> scale_factor = rStartPoint.GetReferenceLengthScale();
376  c_vector<double, DIM> loc0 = locations[0].GetLocation(scale_factor);
377  c_vector<double, DIM> loc1 = locations[1].GetLocation(scale_factor);
378  c_vector<double, DIM> loc2 = locations[2].GetLocation(scale_factor);
379  c_vector<double, DIM> loc3 = locations[3].GetLocation(scale_factor);
380 
381  vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints> :: New();
382  points->InsertNextPoint(&loc0[0]);
383  points->InsertNextPoint(&loc1[0]);
384  points->InsertNextPoint(&loc2[0]);
385  points->InsertNextPoint(&loc3[0]);
386 
387  vtkSmartPointer<vtkUnstructuredGrid> p_grid = vtkSmartPointer<vtkUnstructuredGrid>::New();
388  p_grid->SetPoints(points);
389 
390  vtkIdType ptIds[] = {0, 1, 2, 3};
391  p_grid->InsertNextCell(VTK_TETRA, 4, ptIds);
392 
393  double t;
394  c_vector<double,DIM> intersection;
395  c_vector<double,DIM> parametric_intersection;
396  int subId;
397 
398  c_vector<double,3> dimensionless_start = rStartPoint.GetLocation(scale_factor);
399  c_vector<double,3> dimensionless_end = rEndPoint.GetLocation(scale_factor);
400 
401  if(point1_in_tetra)
402  {
403  p_grid->GetCell(0)->IntersectWithLine(&dimensionless_start[0], &dimensionless_end[0], 1.e-6, t, &intersection[0], &parametric_intersection[0], subId);
404  return norm_2(intersection - dimensionless_start)*scale_factor;
405  }
406 
407  if(point2_in_tetra)
408  {
409  p_grid->GetCell(0)->IntersectWithLine(&dimensionless_end[0], &dimensionless_start[0], 1.e-6, t, &intersection[0], &parametric_intersection[0], subId);
410  return norm_2(intersection - dimensionless_end)*scale_factor;
411  }
412 
413  line_crosses = p_grid->GetCell(0)->IntersectWithLine(&dimensionless_start[0], &dimensionless_end[0], 1.e-6, t, &intersection[0], &parametric_intersection[0], subId);
414  if(line_crosses)
415  {
416  c_vector<double,DIM> intersection2;
417  p_grid->GetCell(0)->IntersectWithLine(&dimensionless_end[0], &dimensionless_start[0], 1.e-6, t, &intersection2[0], &parametric_intersection[0], subId);
418  return norm_2(intersection - intersection2)*scale_factor;
419  }
420  else
421  {
422  return 0.0*scale_factor;
423  }
424  }
425 }
426 
427 template<unsigned DIM>
428 DimensionalChastePoint<DIM> OffsetAlongVector(const DimensionalChastePoint<DIM>& rVector, units::quantity<unit::length> offset)
429 {
430  return DimensionalChastePoint<DIM>(rVector.GetUnitVector() * double(offset/rVector.GetReferenceLengthScale()), rVector.GetReferenceLengthScale());
431 }
432 
433 template<unsigned DIM>
434 DimensionalChastePoint<DIM> OffsetAlongVector(const c_vector<double, DIM>& rVector, units::quantity<unit::length> offset,
435  units::quantity<unit::length> referenceLength)
436 {
437  return DimensionalChastePoint<DIM>(rVector * double(offset/referenceLength), referenceLength);
438 }
439 
440 template<unsigned DIM>
441 DimensionalChastePoint<DIM> RotateAboutAxis(const DimensionalChastePoint<DIM>& rDirection,
442  const DimensionalChastePoint<3>& rAxis, units::quantity<unit::plane_angle> angle)
443 {
444  double sin_a = units::sin(angle);
445  double cos_a = units::cos(angle);
446  c_vector<double, DIM> new_direction;
447  units::quantity<unit::length> length_scale = rDirection.GetReferenceLengthScale();
448  c_vector<double, DIM> dimensionless_direction = rDirection.GetLocation(length_scale);
449  if(DIM==3)
450  {
451  c_vector<double, 3> unit_axis = rAxis.GetUnitVector();
452  units::quantity<unit::length> dot_product = GetDotProduct(rDirection, unit_axis);
453  double dimensionless_dot_product = dot_product/length_scale;
454  new_direction[0] = (unit_axis[0] * dimensionless_dot_product * (1.0 - cos_a) + dimensionless_direction[0] * cos_a
455  + (-unit_axis[2] * dimensionless_direction[1] + unit_axis[1] * dimensionless_direction[2]) * sin_a);
456  new_direction[1] = (unit_axis[1] * dimensionless_dot_product * (1.0 - cos_a) + dimensionless_direction[1] * cos_a
457  + (unit_axis[2] * dimensionless_direction[0] - unit_axis[0] * dimensionless_direction[2]) * sin_a);
458  new_direction[2] = (unit_axis[2] * dimensionless_dot_product * (1.0 - cos_a) + dimensionless_direction[2] * cos_a
459  + (-unit_axis[1] * dimensionless_direction[0] + unit_axis[0] * dimensionless_direction[1]) * sin_a);
460  }
461  else
462  {
463  new_direction[0] = dimensionless_direction[0] * cos_a - dimensionless_direction[1]*sin_a;
464  new_direction[1] = dimensionless_direction[0] * sin_a + dimensionless_direction[1]*cos_a;
465  }
466  return DimensionalChastePoint<DIM>(new_direction, length_scale);
467 }
468 
469 template<unsigned DIM>
470 c_vector<double, DIM> RotateAboutAxis(const c_vector<double, DIM>& rDirection,
471  const c_vector<double, 3>& rAxis, units::quantity<unit::plane_angle> angle)
472 {
473  double sin_a = units::sin(angle);
474  double cos_a = units::cos(angle);
475  c_vector<double, DIM> new_direction = zero_vector<double>(DIM);
476 
477  if(DIM==3)
478  {
479  c_vector<double, 3> unit_axis = rAxis/norm_2(rAxis);
480  double dot_product = inner_prod(rDirection, unit_axis);
481  new_direction[0] = (unit_axis[0] * dot_product * (1.0 - cos_a) + rDirection[0] * cos_a
482  + (-unit_axis[2] * rDirection[1] + unit_axis[1] * rDirection[2]) * sin_a);
483  new_direction[1] = (unit_axis[1] * dot_product * (1.0 - cos_a) + rDirection[1] * cos_a
484  + (unit_axis[2] * rDirection[0] - unit_axis[0] * rDirection[2]) * sin_a);
485  new_direction[2] = (unit_axis[2] * dot_product * (1.0 - cos_a) + rDirection[2] * cos_a
486  + (-unit_axis[1] * rDirection[0] + unit_axis[0] * rDirection[1]) * sin_a);
487  }
488  else
489  {
490  new_direction[0] = rDirection[0] * cos_a - rDirection[1]*sin_a;
491  new_direction[1] = rDirection[0] * sin_a + rDirection[1]*cos_a;
492  }
493  return new_direction;
494 }
units::quantity< unit::length > GetReferenceLengthScale() const
units::quantity< unit::length > GetDistance(const DimensionalChastePoint< DIM > &rLocation) const
c_vector< double, DIM > GetUnitVector() const
c_vector< double, DIM > GetLocation(units::quantity< unit::length > scale)
units::quantity< unit::length > GetNorm2()