Chaste  Build::
DimensionalChastePoint.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 "UblasIncludes.hpp"
37 #include "DimensionalChastePoint.hpp"
38 
39 template<unsigned DIM>
40 DimensionalChastePoint<DIM>::DimensionalChastePoint(double x, double y, double z, units::quantity<unit::length> referenceLength) :
41  mReferenceLength(referenceLength),
42  mIndex(0)
43 {
44  if(mReferenceLength == 0.0*unit::metres)
45  {
46  EXCEPTION("Point has zero reference length");
47  }
48 
49  if (DIM > 0)
50  {
51  mLocation[0] = x;
52  }
53  if (DIM > 1)
54  {
55  mLocation[1] = y;
56  }
57  if (DIM > 2)
58  {
59  mLocation[2] = z;
60  }
61 }
62 
63 template<unsigned DIM>
64 DimensionalChastePoint<DIM>::DimensionalChastePoint(c_vector<double, DIM> coords, units::quantity<unit::length> referenceLength) :
65  mReferenceLength(referenceLength),
66  mIndex(0)
67 {
68  if(mReferenceLength == 0.0*unit::metres)
69  {
70  EXCEPTION("Point has zero reference length");
71  }
72 
73  for (unsigned i=0; i<DIM; i++)
74  {
75  mLocation(i) = coords[i];
76  }
77 }
78 
79 template<unsigned DIM>
80 boost::shared_ptr<DimensionalChastePoint<DIM> > DimensionalChastePoint<DIM>::Create(double x, double y, double z, units::quantity<unit::length> referenceLength)
81 {
82  MAKE_PTR_ARGS(DimensionalChastePoint<DIM>, p_point, (x, y, z, referenceLength));
83  return p_point;
84 }
85 
86 template<unsigned DIM>
87 boost::shared_ptr<DimensionalChastePoint<DIM> > DimensionalChastePoint<DIM>::Create(c_vector<double, DIM> coords, units::quantity<unit::length> referenceLength)
88 {
89  MAKE_PTR_ARGS(DimensionalChastePoint<DIM>, p_point, (coords, referenceLength));
90  return p_point;
91 }
92 
93 template<unsigned DIM>
95 {
96 
97 }
98 
99 template<unsigned DIM>
100 c_vector<double, DIM> DimensionalChastePoint<DIM>::GetLocation(units::quantity<unit::length> scale)
101 {
102  return mLocation*(mReferenceLength/scale);
103 }
104 
105 template<unsigned DIM>
106 const c_vector<double, DIM> DimensionalChastePoint<DIM>::GetLocation(units::quantity<unit::length> scale) const
107 {
108  return mLocation*(mReferenceLength/scale);
109 }
110 
111 template<unsigned DIM>
112 units::quantity<unit::length> DimensionalChastePoint<DIM>::GetReferenceLengthScale() const
113 {
114  return mReferenceLength;
115 }
116 
117 template<unsigned DIM>
118 units::quantity<unit::length> DimensionalChastePoint<DIM>::GetDistance(const DimensionalChastePoint<DIM>& rLocation) const
119 {
120  return norm_2(rLocation.GetLocation(mReferenceLength) - mLocation)*mReferenceLength;
121 }
122 
123 template<unsigned DIM>
124 units::quantity<unit::length> DimensionalChastePoint<DIM>::GetNorm2()
125 {
126  return norm_2(mLocation)*mReferenceLength;
127 }
128 
129 template<unsigned DIM>
130 c_vector<double, DIM> DimensionalChastePoint<DIM>::GetUnitVector() const
131 {
132  return mLocation/norm_2(mLocation);
133 }
134 
135 template<unsigned DIM>
137 {
139 }
140 
141 template<unsigned DIM>
143 {
144  return (rLocation.GetLocation(mReferenceLength)- mLocation) / (GetDistance(rLocation)/mReferenceLength);
145 }
146 
147 template<unsigned DIM>
149 {
150  mLocation /= factor;
151  return *this;
152 }
153 
154 template<unsigned DIM>
156 {
157  mLocation *= factor;
158  return *this;
159 }
160 
161 template<unsigned DIM>
163 {
164  mLocation += rLocation.GetLocation(mReferenceLength);
165  return *this;
166 }
167 
168 template<unsigned DIM>
170 {
171  mLocation -= rLocation.GetLocation(mReferenceLength);
172  return *this;
173 }
174 
175 template<unsigned DIM>
177 {
178  bool returned_value = true;
179  c_vector<double, DIM> comparison_loc = rLocation.GetLocation(mReferenceLength);
180  for (unsigned dim=0; dim<DIM; dim++)
181  {
182  if (comparison_loc[dim] != mLocation[dim])
183  {
184  returned_value = false;
185  break;
186  }
187  }
188  return returned_value;
189 }
190 
191 template<unsigned DIM>
192 void DimensionalChastePoint<DIM>::SetReferenceLengthScale(units::quantity<unit::length> lenthScale)
193 {
194  if(lenthScale == 0.0*unit::metres)
195  {
196  EXCEPTION("Attempted to assign a zero length scale");
197  }
198  mLocation *= (mReferenceLength/lenthScale);
199  mReferenceLength = lenthScale;
200 }
201 
202 template<unsigned DIM>
203 void DimensionalChastePoint<DIM>::RotateAboutAxis(c_vector<double, 3> axis, double angle)
204 {
205  double sin_a = std::sin(angle);
206  double cos_a = std::cos(angle);
207  c_vector<double, 3> unit_axis = axis / norm_2(axis);
208  if(DIM==2 and unit_axis[2]!= 1.0)
209  {
210  EXCEPTION("2D rotation is about z axis only");
211  }
212 
213  c_vector<double, DIM> old_location = this->mLocation;
214  c_vector<double, DIM> new_location;
215  if(DIM==3)
216  {
217  double dot_product = inner_prod(old_location, unit_axis);
218  new_location[0] = (unit_axis[0] * dot_product * (1.0 - cos_a) + old_location[0] * cos_a
219  + (-unit_axis[2] * old_location[1] + unit_axis[1] * old_location[2]) * sin_a);
220  new_location[1] = (unit_axis[1] * dot_product * (1.0 - cos_a) + old_location[1] * cos_a
221  + (unit_axis[2] * old_location[0] - unit_axis[0] * old_location[2]) * sin_a);
222  new_location[2] = (unit_axis[2] * dot_product * (1.0 - cos_a) + old_location[2] * cos_a
223  + (-unit_axis[1] * old_location[0] + unit_axis[0] * old_location[1]) * sin_a);
224  }
225  else
226  {
227  new_location[0] = old_location[0] * cos_a - old_location[1]*sin_a;
228  new_location[1] = old_location[0] * sin_a + old_location[1]*cos_a;
229  }
230 
231  this->mLocation = new_location;
232 }
233 
234 template<unsigned DIM>
236 {
238 }
239 
240 template<unsigned DIM>
242 {
244 }
245 
246 template<unsigned DIM>
248 {
249  return mIndex;
250 }
251 
252 template<unsigned DIM>
254 {
255  mIndex = index;
256 }
257 
258 // Explicit instantiation
259 template class DimensionalChastePoint<2>;
260 template class DimensionalChastePoint<3>;
DimensionalChastePoint(double x=0.0, double y=0.0, double z=0.0, units::quantity< unit::length > referenceLength=1.e-6 *unit::metres)
c_vector< double, DIM > GetUnitTangent(const DimensionalChastePoint< DIM > &rLocation) const
DimensionalChastePoint< DIM > & operator-=(const DimensionalChastePoint< DIM > &rLocation)
bool IsCoincident(const DimensionalChastePoint< DIM > &rLocation) const
units::quantity< unit::length > GetReferenceLengthScale() const
DimensionalChastePoint< DIM > GetMidPoint(const DimensionalChastePoint< DIM > &rLocation) const
units::quantity< unit::length > GetDistance(const DimensionalChastePoint< DIM > &rLocation) const
DimensionalChastePoint< DIM > & operator*=(double factor)
void TranslateTo(DimensionalChastePoint< DIM > rPoint)
c_vector< double, DIM > GetUnitVector() const
c_vector< double, DIM > GetLocation(units::quantity< unit::length > scale)
units::quantity< unit::length > mReferenceLength
void RotateAboutAxis(c_vector< double, 3 > axis, double angle)
units::quantity< unit::length > GetNorm2()
void SetReferenceLengthScale(units::quantity< unit::length > lenthScale)
static boost::shared_ptr< DimensionalChastePoint< DIM > > Create(double x, double y, double z, units::quantity< unit::length > referenceLength)
DimensionalChastePoint< DIM > & operator/=(double factor)
void Translate(DimensionalChastePoint< DIM > rVector)
c_vector< double, DIM > mLocation
DimensionalChastePoint< DIM > & operator+=(const DimensionalChastePoint< DIM > &rLocation)