Transferred copyright to the OpenFOAM Foundation
[OpenFOAM-2.0.x.git] / src / OpenFOAM / fields / Fields / transformField / transformField.C
blob362ef156ca105b9a189a263bb5a8cff6b1aa76d4
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.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "transformField.H"
27 #include "FieldM.H"
28 #include "diagTensor.H"
30 // * * * * * * * * * * * * * * * global functions  * * * * * * * * * * * * * //
32 void Foam::transform
34     vectorField& rtf,
35     const quaternion& q,
36     const vectorField& tf
39     tensor t = q.R();
40     TFOR_ALL_F_OP_FUNC_S_F(vector, rtf, =, transform, tensor, t, vector, tf)
44 Foam::tmp<Foam::vectorField> Foam::transform
46     const quaternion& q,
47     const vectorField& tf
50     tmp<vectorField > tranf(new vectorField(tf.size()));
51     transform(tranf(), q, tf);
52     return tranf;
56 Foam::tmp<Foam::vectorField> Foam::transform
58     const quaternion& q,
59     const tmp<vectorField>& ttf
62     tmp<vectorField > tranf = reuseTmp<vector, vector>::New(ttf);
63     transform(tranf(), q, ttf());
64     reuseTmp<vector, vector>::clear(ttf);
65     return tranf;
69 void Foam::transform
71     vectorField& rtf,
72     const septernion& tr,
73     const vectorField& tf
76     vector T = tr.t();
78     // Check if any rotation
79     if (mag(tr.r().R() - I) > SMALL)
80     {
81         transform(rtf, tr.r(), tf);
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     }
102 Foam::tmp<Foam::vectorField> Foam::transform
104     const septernion& tr,
105     const vectorField& tf
108     tmp<vectorField > tranf(new vectorField(tf.size()));
109     transform(tranf(), tr, tf);
110     return tranf;
114 Foam::tmp<Foam::vectorField> Foam::transform
116     const septernion& tr,
117     const tmp<vectorField>& ttf
120     tmp<vectorField > tranf = reuseTmp<vector, vector>::New(ttf);
121     transform(tranf(), tr, ttf());
122     reuseTmp<vector, vector>::clear(ttf);
123     return tranf;
127 // ************************************************************************* //