motionSmootherAlgoTemplates.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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "motionSmootherAlgo.H"
27 #include "meshTools.H"
29 #include "pointConstraint.H"
30 #include "pointConstraints.H"
31 #include "syncTools.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
39 )
40 {
42 
43  const polyMesh& mesh = pf.mesh();
44 
45  const polyBoundaryMesh& bm = mesh.boundaryMesh();
46 
47  // first count the total number of patch-patch points
48 
49  label nPatchPatchPoints = 0;
50 
51  forAll(bm, patchi)
52  {
53  if (!isA<emptyPolyPatch>(bm[patchi]))
54  {
55  nPatchPatchPoints += bm[patchi].boundaryPoints().size();
56  }
57  }
58 
59 
60  typename FldType::GeometricBoundaryField& bFld = pf.boundaryField();
61 
62 
63  // Evaluate in reverse order
64 
65  forAllReverse(bFld, patchi)
66  {
67  bFld[patchi].initEvaluate(Pstream::blocking); // buffered
68  }
69 
70  forAllReverse(bFld, patchi)
71  {
72  bFld[patchi].evaluate(Pstream::blocking);
73  }
74 
75 
76  // Save the values
77 
78  Field<Type> boundaryPointValues(nPatchPatchPoints);
79  nPatchPatchPoints = 0;
80 
81  forAll(bm, patchi)
82  {
83  if (!isA<emptyPolyPatch>(bm[patchi]))
84  {
85  const labelList& bp = bm[patchi].boundaryPoints();
86  const labelList& meshPoints = bm[patchi].meshPoints();
87 
88  forAll(bp, pointi)
89  {
90  label ppp = meshPoints[bp[pointi]];
91  boundaryPointValues[nPatchPatchPoints++] = pf[ppp];
92  }
93  }
94  }
95 
96 
97  // Forward evaluation
98 
99  bFld.evaluate();
100 
101 
102  // Check
103 
104  nPatchPatchPoints = 0;
105 
106  forAll(bm, patchi)
107  {
108  if (!isA<emptyPolyPatch>(bm[patchi]))
109  {
110  const labelList& bp = bm[patchi].boundaryPoints();
111  const labelList& meshPoints = bm[patchi].meshPoints();
112 
113  forAll(bp, pointi)
114  {
115  label ppp = meshPoints[bp[pointi]];
116 
117  const Type& savedVal = boundaryPointValues[nPatchPatchPoints++];
118 
119  if (savedVal != pf[ppp])
120  {
122  << "Patch fields are not consistent on mesh point "
123  << ppp << " coordinate " << mesh.points()[ppp]
124  << " at patch " << bm[patchi].name() << '.'
125  << endl
126  << "Reverse evaluation gives value " << savedVal
127  << " , forward evaluation gives value " << pf[ppp]
128  << abort(FatalError);
129  }
130  }
131  }
132  }
133 }
134 
135 
136 // Average of connected points.
137 template<class Type>
140 (
142  const scalarField& edgeWeight
143 ) const
144 {
146  (
148  (
149  IOobject
150  (
151  "avg("+fld.name()+')',
152  fld.time().timeName(),
153  fld.db(),
154  IOobject::NO_READ,
155  IOobject::NO_WRITE,
156  false
157  ),
158  fld.mesh(),
159  dimensioned<Type>("zero", fld.dimensions(), pTraits<Type>::zero)
160  )
161  );
163 
164  const polyMesh& mesh = fld.mesh()();
165 
166 
167  // Sum local weighted values and weights
168  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
169 
170  // Note: on coupled edges use only one edge (through isMasterEdge)
171  // This is done so coupled edges do not get counted double.
172 
173  scalarField sumWeight(mesh.nPoints(), 0.0);
174 
175  const edgeList& edges = mesh.edges();
176 
177  forAll(edges, edgeI)
178  {
179  if (isMasterEdge_.get(edgeI) == 1)
180  {
181  const edge& e = edges[edgeI];
182  const scalar w = edgeWeight[edgeI];
183 
184  res[e[0]] += w*fld[e[1]];
185  sumWeight[e[0]] += w;
186 
187  res[e[1]] += w*fld[e[0]];
188  sumWeight[e[1]] += w;
189  }
190  }
191 
192 
193  // Add coupled contributions
194  // ~~~~~~~~~~~~~~~~~~~~~~~~~
195 
196  syncTools::syncPointList
197  (
198  mesh,
199  res,
200  plusEqOp<Type>(),
201  pTraits<Type>::zero // null value
202  );
203  syncTools::syncPointList
204  (
205  mesh,
206  sumWeight,
208  scalar(0) // null value
209  );
210 
211 
212  // Average
213  // ~~~~~~~
214 
215  forAll(res, pointI)
216  {
217  if (mag(sumWeight[pointI]) < VSMALL)
218  {
219  // Unconnected point. Take over original value
220  res[pointI] = fld[pointI];
221  }
222  else
223  {
224  res[pointI] /= sumWeight[pointI];
225  }
226  }
227 
228  // Single and multi-patch constraints
229  pointConstraints::New(fld.mesh()).constrain(res, false);
230 
231  return tres;
232 }
233 
234 
235 // smooth field (point-jacobi)
236 template<class Type>
238 (
240  const scalarField& edgeWeight,
242 ) const
243 {
244  tmp<pointVectorField> tavgFld = avg(fld, edgeWeight);
245  const pointVectorField& avgFld = tavgFld();
246 
247  forAll(fld, pointI)
248  {
249  if (isInternalPoint(pointI))
250  {
251  newFld[pointI] = 0.5*fld[pointI] + 0.5*avgFld[pointI];
252  }
253  }
254 
255  // Single and multi-patch constraints
256  pointConstraints::New(fld.mesh()).constrain(newFld, false);
257 }
258 
259 
260 //- Test synchronisation of generic field (not positions!) on points
261 template<class Type, class CombineOp>
263 (
264  const Field<Type>& fld,
265  const CombineOp& cop,
266  const Type& zero,
267  const scalar maxMag
268 ) const
269 {
270  if (debug)
271  {
272  Pout<< "testSyncField : testing synchronisation of Field<Type>."
273  << endl;
274  }
275 
276  Field<Type> syncedFld(fld);
277 
278  syncTools::syncPointList
279  (
280  mesh_,
281  syncedFld,
282  cop,
283  zero
284  );
285 
286  forAll(syncedFld, i)
287  {
288  if (mag(syncedFld[i] - fld[i]) > maxMag)
289  {
291  << "On element " << i << " value:" << fld[i]
292  << " synchronised value:" << syncedFld[i]
293  << abort(FatalError);
294  }
295  }
296 }
297 
298 
299 // ************************************************************************* //
pointConstraint.H
Foam::motionSmootherAlgo::avg
tmp< GeometricField< Type, pointPatchField, pointMesh > > avg(const GeometricField< Type, pointPatchField, pointMesh > &fld, const scalarField &edgeWeight) const
Average of connected points.
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
meshTools.H
w
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::motionSmootherAlgo::testSyncField
void testSyncField(const Field< Type > &, const CombineOp &cop, const Type &zero, const scalar maxMag) const
Test synchronisation of generic field (not positions!) on points.
Definition: motionSmootherAlgoTemplates.C:263
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field< Type >
pointConstraints.H
Foam::motionSmootherAlgo::checkConstraints
static void checkConstraints(GeometricField< Type, pointPatchField, pointMesh > &)
Check constraints.
Definition: motionSmootherAlgoTemplates.C:37
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::FatalError
error FatalError
forAllReverse
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:418
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::plusEqOp
Definition: ops.H:71
Foam::dimensioned< Type >
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
processorPointPatchFields.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
motionSmootherAlgo.H
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::motionSmootherAlgo::smooth
void smooth(const GeometricField< Type, pointPatchField, pointMesh > &fld, const scalarField &edgeWeight, GeometricField< Type, pointPatchField, pointMesh > &newFld) const
Fully explicit smoothing of fields (not positions)
Definition: motionSmootherAlgoTemplates.C:238
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::zero
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:47