All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
transformField.C
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 #include "transformField.H"
27 #include "FieldM.H"
28 #include "diagTensor.H"
29 
30 // * * * * * * * * * * * * * * * global functions * * * * * * * * * * * * * //
31 
32 void Foam::transform
33 (
34  vectorField& rtf,
35  const quaternion& q,
36  const vectorField& tf
37 )
38 {
39  tensor t = q.R();
41 }
42 
43 
45 (
46  const quaternion& q,
47  const vectorField& tf
48 )
49 {
50  tmp<vectorField > tranf(new vectorField(tf.size()));
51  transform(tranf.ref(), q, tf);
52  return tranf;
53 }
54 
55 
57 (
58  const quaternion& q,
59  const tmp<vectorField>& ttf
60 )
61 {
62  tmp<vectorField > tranf = New(ttf);
63  transform(tranf.ref(), q, ttf());
64  ttf.clear();
65  return tranf;
66 }
67 
68 
70 (
71  vectorField& rtf,
72  const septernion& tr,
73  const vectorField& tf
74 )
75 {
76  vector T = tr.t();
77 
78  // Check if any translation
79  if (mag(T) > vSmall)
80  {
81  TFOR_ALL_F_OP_F_OP_S(vector, rtf, =, vector, tf, -, vector, T);
82  }
83  else
84  {
85  rtf = tf;
86  }
87 
88  // Check if any rotation
89  if (mag(tr.r().R() - I) > small)
90  {
91  transform(rtf, tr.r(), rtf);
92  }
93 }
94 
95 
97 (
98  const septernion& tr,
99  const vectorField& tf
100 )
101 {
102  tmp<vectorField > tranf(new vectorField(tf.size()));
103  transformPoints(tranf.ref(), tr, tf);
104  return tranf;
105 }
106 
107 
109 (
110  const septernion& tr,
111  const tmp<vectorField>& ttf
112 )
113 {
114  tmp<vectorField > tranf = New(ttf);
115  transformPoints(tranf.ref(), tr, ttf());
116  ttf.clear();
117  return tranf;
118 }
119 
120 
121 // ************************************************************************* //
tensor R() const
The rotation tensor corresponding the quaternion.
Definition: quaternionI.H:323
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
void transformPoints(vectorField &, const septernion &, const vectorField &)
Transform given vectorField of coordinates with the given septernion.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
volVectorField vectorField(fieldObject, mesh)
Septernion class used to perform translations and rotations in 3D space.
Definition: septernion.H:65
const tensorField & tf
Spatial transformation functions for primitive fields.
static const Identity< scalar > I
Definition: Identity.H:93
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:60
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
const volScalarField & T
#define TFOR_ALL_F_OP_FUNC_S_F(typeF1, f1, OP, FUNC, typeS, s, typeF2, f2)
Definition: FieldM.H:210
const vector & t() const
Definition: septernionI.H:58
dimensioned< scalar > mag(const dimensioned< Type > &)
const quaternion & r() const
Definition: septernionI.H:64
A class for managing temporary objects.
Definition: PtrList.H:53
High performance macro functions for Field<Type> algebra. These expand using either array element acc...
#define TFOR_ALL_F_OP_F_OP_S(typeF1, f1, OP1, typeF2, f2, OP2, typeS, s)
Definition: FieldM.H:295
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477