transformField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011 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(), q, tf);
52  return tranf;
53 }
54 
55 
57 (
58  const quaternion& q,
59  const tmp<vectorField>& ttf
60 )
61 {
63  transform(tranf(), q, ttf());
65  return tranf;
66 }
67 
68 
69 void Foam::transform
70 (
71  vectorField& rtf,
72  const septernion& tr,
73  const vectorField& tf
74 )
75 {
76  vector T = tr.t();
77 
78  // Check if any rotation
79  if (mag(tr.r().R() - I) > SMALL)
80  {
81  transform(rtf, tr.r(), tf);
82 
83  if (mag(T) > VSMALL)
84  {
85  rtf += T;
86  }
87  }
88  else
89  {
90  if (mag(T) > VSMALL)
91  {
92  TFOR_ALL_F_OP_S_OP_F(vector, rtf, =, vector, T, +, vector, tf);
93  }
94  else
95  {
96  rtf = tf;
97  }
98  }
99 }
100 
101 
103 (
104  const septernion& tr,
105  const vectorField& tf
106 )
107 {
108  tmp<vectorField > tranf(new vectorField(tf.size()));
109  transform(tranf(), tr, tf);
110  return tranf;
111 }
112 
113 
115 (
116  const septernion& tr,
117  const tmp<vectorField>& ttf
118 )
119 {
121  transform(tranf(), tr, ttf());
123  return tranf;
124 }
125 
126 
127 // ************************************************************************* //
Foam::Tensor
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
Definition: complexI.H:224
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::septernion
Septernion class used to perform translations and rotations in 3D space.
Definition: septernion.H:64
clear
UEqn clear()
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
FieldM.H
High performance macro functions for Field<Type> algebra. These expand using either array element acc...
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
transformField.H
Spatial transformation functions for primitive fields.
Foam::transform
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
Foam::quaternion
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:60
tf
const tensorField & tf
Definition: getPatchFieldTensor.H:36
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
TFOR_ALL_F_OP_S_OP_F
#define TFOR_ALL_F_OP_S_OP_F(typeF1, f1, OP1, typeS, s, OP2, typeF2, f2)
Definition: FieldM.H:278
vectorField
volVectorField vectorField(fieldObject, mesh)
Foam::I
static const sphericalTensor I(1)
T
const volScalarField & T
Definition: createFields.H:25
TFOR_ALL_F_OP_FUNC_S_F
#define TFOR_ALL_F_OP_FUNC_S_F(typeF1, f1, OP, FUNC, typeS, s, typeF2, f2)
Definition: FieldM.H:210
Foam::Vector< scalar >
Foam::quaternion::R
tensor R() const
The rotation tensor corresponding the quaternion.
Definition: quaternionI.H:235
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:49
diagTensor.H