vectorField.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "vectorField.H"
29
30// * * * * * * * * * * * * * * * Specializations * * * * * * * * * * * * * * //
31
32namespace Foam
33{
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37// Note: the mag of vector and division is written out to avoid any casting
38// between float and double.
39//
40// This enables specialization for floatVector and doubleVector independent
41// of the definition of 'scalar' or 'vector' - useful for mixed-precision
42// operation.
43
44template<>
45void Field<Vector<float>>::normalise()
46{
47 typedef float cmptType;
48
49 constexpr float tol = floatScalarROOTVSMALL;
50
51 for (Vector<cmptType>& v : *this)
52 {
53 // Foam::mag but using cmptType instead of scalar
55 (
56 ::sqrt(magSqr(v.x()) + magSqr(v.y()) + magSqr(v.z()))
57 );
58
59 if (s < tol)
60 {
61 v.x() = 0;
62 v.y() = 0;
63 v.z() = 0;
64 }
65 else
66 {
67 v.x() /= s;
68 v.y() /= s;
69 v.z() /= s;
70 }
71 }
72}
73
74
75template<>
76void Field<Vector<double>>::normalise()
77{
78 typedef double cmptType;
79
80 constexpr double tol = doubleScalarROOTVSMALL;
81
82 for (Vector<cmptType>& v : *this)
83 {
84 // Foam::mag but using cmptType instead of scalar
86 (
87 ::sqrt(magSqr(v.x()) + magSqr(v.y()) + magSqr(v.z()))
88 );
89
90 if (s < tol)
91 {
92 v.x() = 0;
93 v.y() = 0;
94 v.z() = 0;
95 }
96 else
97 {
98 v.x() /= s;
99 v.y() /= s;
100 v.z() /= s;
101 }
102 }
103}
104
105
106// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
107
108} // End namespace Foam
109
110// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
111
112// ************************************************************************* //
Generic templated field type.
Definition: Field.H:82
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:65
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
constexpr floatScalar floatScalarROOTVSMALL
Definition: floatScalar.H:64
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr doubleScalar doubleScalarROOTVSMALL
Definition: doubleScalar.H:64