weightedPosition.H
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
26Class
27 Foam::weightedPosition
28
29Description
30 Wrapper for position + weight to be used in e.g. averaging.
31
32 This avoids the problems when synchronising locations with e.g.
33 parallel cyclics. The separation vector only applies to locations
34 and not e.g. summed locations. Note that there is no corresponding
35 problem for rotational cyclics.
36
37 Typical use might be to e.g. average face centres to points on a patch
38
39 const labelListList& pointFaces = pp.pointFaces();
40 const vectorField::subField faceCentres = pp.faceCentres();
41
42 Field<weightedPosition> avgBoundary(pointFaces.size());
43
44 forAll(pointFaces, pointi)
45 {
46 const labelList& pFaces = pointFaces[pointi];
47 avgBoundary[pointi].first() = pFaces.size();
48 avgBoundary[pointi].second() = sum(pointField(faceCentres, pFaces));
49 }
50 syncTools::syncPointList
51 (
52 mesh,
53 pp.meshPoints(),
54 avgBoundary,
55 weightedPosition::plusEqOp, // combine op
56 pTraits<weightedPosition>::zero,// null value (not used)
57 pTraits<weightedPosition>::zero // transform class
58 );
59
60
61SourceFiles
62 weightedPosition.C
63
64\*---------------------------------------------------------------------------*/
65
66#ifndef weightedPosition_H
67#define weightedPosition_H
68
69#include "Tuple2.H"
70#include "point.H"
71#include "Field.H"
72
73// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74
75namespace Foam
76{
77
78// Forward declarations
79class weightedPosition;
80class vectorTensorTransform;
81class coupledPolyPatch;
82class polyMesh;
83
84//- pTraits
85template<>
87{
88public:
90 static const weightedPosition zero;
91};
92
93
94/*---------------------------------------------------------------------------*\
95 Class weightedPosition Declaration
96\*---------------------------------------------------------------------------*/
99:
100 public Tuple2<scalar, point>
101{
102public:
103
104 // Constructors
105
106 //- Construct null
108
109 //- Construct from components
110 weightedPosition(const scalar s, const point& p);
111
112
113 // Member Functions
114
115 // Helpers
116
117 //- Get points
118 static void getPoints
119 (
120 const UList<weightedPosition>& in,
121 List<point>& out
122 );
123
124 //- Set points
125 static void setPoints
126 (
127 const UList<point>& in,
129 );
130
131
132 //- Summation operator
133 static void plusEqOp(weightedPosition& x, const weightedPosition& y);
134
135
136 // Transformation operators
137
138 void operator()
139 (
140 const vectorTensorTransform& vt,
141 const bool forward,
143 ) const;
144
145 void operator()
146 (
147 const vectorTensorTransform& vt,
148 const bool forward,
150 ) const;
151
152 void operator()
153 (
154 const coupledPolyPatch& cpp,
156 ) const;
157
158 template<template<class> class Container>
159 void operator()
160 (
161 const coupledPolyPatch& cpp,
162 Container<weightedPosition>& map
163 ) const;
164
165
166 //- Synchronisation for mesh point positions
167 static void syncPoints
168 (
169 const polyMesh& mesh,
171 );
172
173 //- Synchronisation for patch point positions
174 static void syncPoints
175 (
176 const polyMesh& mesh,
177 const labelUList& meshPoints,
179 );
180};
181
182
183// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184
185} // End namespace Foam
186
187// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188
189#ifdef NoRepository
191#endif
192
193// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194
195#endif
196
197// ************************************************************************* //
scalar y
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic templated field type.
Definition: Field.H:82
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: Tuple2.H:58
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:94
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.
Definition: pTraits.H:59
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
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.
weightedPosition()
Construct null.
static void syncPoints(const polyMesh &mesh, List< weightedPosition > &)
Synchronisation for mesh point positions.
static void setPoints(const UList< point > &in, List< weightedPosition > &out)
Set points.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
volScalarField & p
dynamicFvMesh & mesh
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.