pointConstraintI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
29 :
30  Tuple2<label, vector>(0, Zero)
31 {}
32 
33 
35 :
36  Tuple2<label, vector>(pc)
37 {}
38 
39 
41 :
42  Tuple2<label, vector>(is)
43 {}
44 
45 
46 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
47 
49 {
50  if (first() == 0)
51  {
52  first() = 1;
53  second() = cd;
54  }
55  else if (first() == 1)
56  {
57  vector planeNormal = cd ^ second();
58  scalar magPlaneNormal = mag(planeNormal);
59 
60  if (magPlaneNormal > 1e-3)
61  {
62  first() = 2;
63  second() = planeNormal/magPlaneNormal;
64  }
65  }
66  else if (first() == 2)
67  {
68  if (mag(cd & second()) > 1e-3)
69  {
70  first() = 3;
71  second() = Zero;
72  }
73  }
74 }
75 
76 
78 {
79  if (first() == 0)
80  {
81  operator=(pc);
82  }
83  else if (first() == 1)
84  {
85  // Save single normal
86  vector n = second();
87  // Apply to supplied point constaint
88  operator=(pc);
89  applyConstraint(n);
90  }
91  else if (first() == 2)
92  {
93  if (pc.first() == 0)
94  {}
95  else if (pc.first() == 1)
96  {
97  applyConstraint(pc.second());
98  }
99  else if (pc.first() == 2)
100  {
101  // Both constrained to line. Same (+-)direction?
102  if (mag(second() & pc.second()) <= (1.0-1e-3))
103  {
104  // Different directions
105  first() = 3;
106  second() = Zero;
107  }
108  }
109  else
110  {
111  first() = 3;
112  second() = Zero;
113  }
114  }
115 }
116 
117 
119 {
120  if (first() == 0)
121  {
122  return I;
123  }
124  else if (first() == 1)
125  {
126  return I - sqr(second());
127  }
128  else if (first() == 2)
129  {
130  return sqr(second());
131  }
132  else
133  {
134  return Zero;
135  }
136 }
137 
138 
140 (
141  label& n,
142  tensor& tt
143 ) const
144 {
145  n = 3-first();
146 
148 
149  if (first() == 0)
150  {
151  vecs[0] = vector(1, 0, 0);
152  vecs[1] = vector(0, 1, 0);
153  vecs[2] = vector(0, 0, 1);
154  }
155  else if (first() == 1)
156  {
157  const vector& planeDir = second();
158 
159  vecs[0] = vector(1, 0, 0) - planeDir.x()*planeDir;
160 
161  if (mag(vecs[0].x()) < 1e-3)
162  {
163  vecs[0] = vector(0, 1, 0) - planeDir.y()*planeDir;
164  }
165 
166  vecs[0] /= mag(vecs[0]);
167  vecs[1] = vecs[0] ^ planeDir;
168  vecs[1] /= mag(vecs[1]);
169  }
170  else if (first() == 2)
171  {
172  vecs[0] = second();
173  }
174 
175  // Knock out remaining vectors
176  for (direction dir = n; dir < vecs.size(); dir++)
177  {
178  vecs[dir] = Zero;
179  }
180 
181  tt = tensor(vecs[0], vecs[1], vecs[2]);
182 }
183 
184 
186 (
187  const vector& d
188 ) const
189 {
190  vector cd;
191 
192  if (first() == 0)
193  {
194  cd = d;
195  }
196  else if (first() == 1)
197  {
198  // Remove plane normal
199  cd = d-(d&second())*second();
200  }
201  else if (first() == 2)
202  {
203  // Keep line direction only
204  cd = (d&second())*second();
205  }
206  else
207  {
208  cd = Zero;
209  }
210  return cd;
211 }
212 
213 
214 inline void Foam::combineConstraintsEqOp::operator()
215 (
217  const pointConstraint& y
218 ) const
219 {
220  x.combine(y);
221 }
222 
223 
225 (
226  const tensor& tt,
227  const pointConstraint& v
228 )
229 {
230  return pointConstraint
231  (
233  );
234 }
235 
236 
237 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
A 2-tuple for storing two objects of different types.
Definition: HashTable.H:66
const label & first() const
Return first.
Definition: Tuple2.H:94
void unconstrainedDirections(label &n, tensor &vecs) const
Return the accumulated unconstrained directions. Directions.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
uint8_t direction
Definition: direction.H:45
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void combine(const pointConstraint &)
Combine constraints.
const Cmpt & y() const
Definition: VectorI.H:81
void applyConstraint(const vector &cd)
Apply and accumulate the effect of the given constraint direction.
scalar y
static const Identity< scalar > I
Definition: Identity.H:93
static const zero Zero
Definition: zero.H:97
Accumulates point constraints through successive applications of the applyConstraint function...
const Cmpt & x() const
Definition: VectorI.H:75
label size() const
Return the number of elements in the FixedList.
Definition: FixedListI.H:429
tensor constraintTransformation() const
Return the accumulated constraint transformation tensor.
vector constrainDisplacement(const vector &disp) const
Constrain a displacement.
const vector & second() const
Return second.
Definition: Tuple2.H:106
dimensioned< scalar > mag(const dimensioned< Type > &)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
pointConstraint()
Construct null.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477