weightedPosition.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) 2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 "weightedPosition.H"
29 #include "vectorTensorTransform.H"
30 #include "coupledPolyPatch.H"
31 #include "polyMesh.H"
32 #include "syncTools.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
37 (
38  scalar(0),
39  Foam::point::zero
40 );
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 :
47  Tuple2<scalar, point>()
48 {}
49 
50 
52 :
53  Tuple2<scalar, point>(s, p)
54 {}
55 
56 
57 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
58 
60 (
61  const UList<weightedPosition>& in,
62  List<point>& out
63 )
64 {
65  out.setSize(in.size());
66  forAll(in, i)
67  {
68  out[i] = in[i].second();
69 
70  if (mag(in[i].first()) > VSMALL)
71  {
72  out[i] /= in[i].first();
73  }
74  }
75 }
76 
77 
79 (
80  const UList<point>& in,
81  List<weightedPosition>& out
82 )
83 {
84  out.setSize(in.size());
85  forAll(in, i)
86  {
87  out[i].second() = out[i].first()*in[i];
88  }
89 }
90 
91 
93 (
94  weightedPosition& x,
95  const weightedPosition& y
96 )
97 {
98  x.first() += y.first();
99  x.second() += y.second();
100 }
101 
102 
103 void Foam::weightedPosition::operator()
104 (
105  const vectorTensorTransform& vt,
106  const bool forward,
108 ) const
109 {
110  pointField pfld;
111  getPoints(fld, pfld);
112 
113  if (forward)
114  {
115  pfld = vt.transformPosition(pfld);
116  }
117  else
118  {
119  pfld = vt.invTransformPosition(pfld);
120  }
121 
122  setPoints(pfld, fld);
123 }
124 
125 
126 void Foam::weightedPosition::operator()
127 (
128  const vectorTensorTransform& vt,
129  const bool forward,
130  List<List<weightedPosition>>& flds
131 ) const
132 {
133  for (List<weightedPosition>& fld : flds)
134  {
135  operator()(vt, forward, fld);
136  }
137 }
138 
139 
140 void Foam::weightedPosition::operator()
141 (
142  const coupledPolyPatch& cpp,
143  Field<weightedPosition>& fld
144 ) const
145 {
146  pointField pfld;
147  getPoints(fld, pfld);
148 
149  cpp.transformPosition(pfld);
150 
151  setPoints(pfld, fld);
152 }
153 
154 
156 (
157  const polyMesh& mesh,
159 )
160 {
161  if (fld.size() != mesh.nPoints())
162  {
163  FatalErrorInFunction << "Size of field " << fld.size()
164  << " does not correspond to the number of points in the mesh "
165  << mesh.nPoints() << exit(FatalError);
166  }
167 
169  (
170  mesh,
171  fld,
172  weightedPosition::plusEqOp, // combine op
173  pTraits<weightedPosition>::zero,// null value (not used)
174  pTraits<weightedPosition>::zero // transform class
175  );
176 }
177 
178 
180 (
181  const polyMesh& mesh,
182  const labelUList& meshPoints,
184 )
185 {
186  if (fld.size() != meshPoints.size())
187  {
188  FatalErrorInFunction << "Size of field " << fld.size()
189  << " does not correspond to the number of points supplied "
190  << meshPoints.size() << exit(FatalError);
191  }
192 
194  (
195  mesh,
196  meshPoints,
197  fld,
198  weightedPosition::plusEqOp, // combine op
199  pTraits<weightedPosition>::zero,// null value (not used)
200  pTraits<weightedPosition>::zero // transform class
201  );
202 }
203 
204 
205 // ************************************************************************* //
weightedPosition()
Construct null.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static void plusEqOp(weightedPosition &x, const weightedPosition &y)
Summation operator.
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:54
static const weightedPosition zero
Vector-tensor class used to perform translations and rotations in 3D space.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
T & first()
Access first element of the list, position [0].
Definition: UList.H:853
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
static void getPoints(const UList< weightedPosition > &in, List< point > &out)
Get points.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
scalar y
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
static void setPoints(const UList< point > &in, List< weightedPosition > &out)
Set points.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Wrapper for position + weight to be used in e.g. averaging.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static void syncPoints(const polyMesh &mesh, List< weightedPosition > &)
Synchronisation for mesh point positions.
volScalarField & p
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))