Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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 "weightedPosition.H"
30#include "coupledPolyPatch.H"
31#include "polyMesh.H"
32#include "syncTools.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
37(
38 scalar(0),
40);
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48{}
49
50
53 Tuple2<scalar, point>(s, p)
54{}
55
56
57// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
58
60(
62 List<point>& out
63)
64{
65 out.resize_nocopy(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,
82)
83{
84 #ifdef FULLDEBUG
85 if (in.size() != out.size())
86 {
88 << "Mismatch in size. Number of points:" << in.size()
89 << " != number of weighted positions:" << out.size()
90 << exit(FatalError);
91 }
92
93 const label len = out.size();
94 #else
95 const label len = Foam::min(in.size(), out.size());
96 #endif
97
98 for (label i = 0; i < len; ++i)
99 {
100 out[i].second() = out[i].first()*in[i];
101 }
102}
103
104
106(
107 weightedPosition& x,
108 const weightedPosition& y
110{
111 x.first() += y.first();
112 x.second() += y.second();
113}
114
115
116void Foam::weightedPosition::operator()
117(
118 const vectorTensorTransform& vt,
119 const bool forward,
121) const
122{
123 pointField pfld;
124 getPoints(fld, pfld);
125
126 if (forward)
127 {
128 vt.transformPositionList(pfld);
129 }
130 else
131 {
132 vt.invTransformPositionList(pfld);
133 }
134
135 setPoints(pfld, fld);
136}
137
138
139void Foam::weightedPosition::operator()
140(
141 const vectorTensorTransform& vt,
142 const bool forward,
143 List<List<weightedPosition>>& flds
144) const
145{
146 for (auto& fld : flds)
147 {
148 operator()(vt, forward, fld);
149 }
150}
151
152
153void Foam::weightedPosition::operator()
154(
155 const coupledPolyPatch& cpp,
156 Field<weightedPosition>& fld
157) const
158{
159 pointField pfld;
160 getPoints(fld, pfld);
162 cpp.transformPosition(pfld);
163
164 setPoints(pfld, fld);
165}
166
167
169(
170 const polyMesh& mesh,
172)
173{
174 if (fld.size() != mesh.nPoints())
175 {
176 FatalErrorInFunction << "Size of field " << fld.size()
177 << " does not correspond to the number of points in the mesh "
178 << mesh.nPoints() << exit(FatalError);
179 }
180
182 (
183 mesh,
184 fld,
186 pTraits<weightedPosition>::zero,// null value (not used)
187 pTraits<weightedPosition>::zero // transform class
188 );
189}
190
191
193(
194 const polyMesh& mesh,
195 const labelUList& meshPoints,
197)
198{
199 if (fld.size() != meshPoints.size())
200 {
201 FatalErrorInFunction << "Size of field " << fld.size()
202 << " does not correspond to the number of points supplied "
203 << meshPoints.size() << exit(FatalError);
204 }
205
207 (
208 mesh,
209 meshPoints,
210 fld,
211 weightedPosition::plusEqOp, // combine op
212 pTraits<weightedPosition>::zero,// null value (not used)
213 pTraits<weightedPosition>::zero // transform class
214 );
215}
216
217
218// ************************************************************************* //
scalar y
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
const scalar & first() const noexcept
Definition Tuple2.H:132
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
T & first()
Access first element of the list, position [0].
Definition UList.H:957
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static const Form zero
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
static const weightedPosition zero
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Vector-tensor class used to perform translations and rotations in 3D space.
Wrapper for position + weight to be used in e.g. averaging.
static void getPoints(const UList< weightedPosition > &in, List< point > &out)
Get points.
static void plusEqOp(weightedPosition &x, const weightedPosition &y)
Summation operator.
weightedPosition()
Construct null.
static void syncPoints(const polyMesh &mesh, List< weightedPosition > &)
Synchronisation for mesh point positions.
static void setPoints(const UList< point > &in, UList< weightedPosition > &out)
Set points.
void operator()(const vectorTensorTransform &vt, const bool forward, UList< weightedPosition > &fld) const
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
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))
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299