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-2013 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>
36 void Foam::motionSmootherAlgo::checkConstraints
37 (
38  GeometricField<Type, pointPatchField, pointMesh>& pf
39 )
40 {
41  typedef GeometricField<Type, pointPatchField, pointMesh> FldType;
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  (
123  "motionSmootherAlgo::checkConstraints"
124  "(GeometricField<Type, pointPatchField, pointMesh>&)"
125  ) << "Patch fields are not consistent on mesh point "
126  << ppp << " coordinate " << mesh.points()[ppp]
127  << " at patch " << bm[patchi].name() << '.'
128  << endl
129  << "Reverse evaluation gives value " << savedVal
130  << " , forward evaluation gives value " << pf[ppp]
131  << abort(FatalError);
132  }
133  }
134  }
135  }
136 }
137 
138 
139 // Average of connected points.
140 template<class Type>
142 Foam::motionSmootherAlgo::avg
143 (
145  const scalarField& edgeWeight
146 ) const
147 {
149  (
151  (
152  IOobject
153  (
154  "avg("+fld.name()+')',
155  fld.time().timeName(),
156  fld.db(),
159  false
160  ),
161  fld.mesh(),
163  )
164  );
166 
167  const polyMesh& mesh = fld.mesh()();
168 
169 
170  // Sum local weighted values and weights
171  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
172 
173  // Note: on coupled edges use only one edge (through isMasterEdge)
174  // This is done so coupled edges do not get counted double.
175 
176  scalarField sumWeight(mesh.nPoints(), 0.0);
177 
178  const edgeList& edges = mesh.edges();
179 
180  forAll(edges, edgeI)
181  {
182  if (isMasterEdge_.get(edgeI) == 1)
183  {
184  const edge& e = edges[edgeI];
185  const scalar w = edgeWeight[edgeI];
186 
187  res[e[0]] += w*fld[e[1]];
188  sumWeight[e[0]] += w;
189 
190  res[e[1]] += w*fld[e[0]];
191  sumWeight[e[1]] += w;
192  }
193  }
194 
195 
196  // Add coupled contributions
197  // ~~~~~~~~~~~~~~~~~~~~~~~~~
198 
200  (
201  mesh,
202  res,
203  plusEqOp<Type>(),
204  pTraits<Type>::zero // null value
205  );
207  (
208  mesh,
209  sumWeight,
211  scalar(0) // null value
212  );
213 
214 
215  // Average
216  // ~~~~~~~
217 
218  forAll(res, pointI)
219  {
220  if (mag(sumWeight[pointI]) < VSMALL)
221  {
222  // Unconnected point. Take over original value
223  res[pointI] = fld[pointI];
224  }
225  else
226  {
227  res[pointI] /= sumWeight[pointI];
228  }
229  }
230 
231  // Single and multi-patch constraints
232  pointConstraints::New(fld.mesh()).constrain(res, false);
233 
234  return tres;
235 }
236 
237 
238 // smooth field (point-jacobi)
239 template<class Type>
241 (
243  const scalarField& edgeWeight,
245 ) const
246 {
247  tmp<pointVectorField> tavgFld = avg(fld, edgeWeight);
248  const pointVectorField& avgFld = tavgFld();
249 
250  forAll(fld, pointI)
251  {
252  if (isInternalPoint(pointI))
253  {
254  newFld[pointI] = 0.5*fld[pointI] + 0.5*avgFld[pointI];
255  }
256  }
257 
258  // Single and multi-patch constraints
259  pointConstraints::New(fld.mesh()).constrain(newFld, false);
260 }
261 
262 
263 //- Test synchronisation of generic field (not positions!) on points
264 template<class Type, class CombineOp>
265 void Foam::motionSmootherAlgo::testSyncField
266 (
267  const Field<Type>& fld,
268  const CombineOp& cop,
269  const Type& zero,
270  const scalar maxMag
271 ) const
272 {
273  if (debug)
274  {
275  Pout<< "testSyncField : testing synchronisation of Field<Type>."
276  << endl;
277  }
278 
279  Field<Type> syncedFld(fld);
280 
282  (
283  mesh_,
284  syncedFld,
285  cop,
286  zero
287  );
288 
289  forAll(syncedFld, i)
290  {
291  if (mag(syncedFld[i] - fld[i]) > maxMag)
292  {
294  (
295  "motionSmootherAlgo::testSyncField"
296  "(const Field<Type>&, const CombineOp&"
297  ", const Type&, const bool)"
298  ) << "On element " << i << " value:" << fld[i]
299  << " synchronised value:" << syncedFld[i]
300  << abort(FatalError);
301  }
302  }
303 }
304 
305 
306 // ************************************************************************* //
dimensioned< scalar > mag(const dimensioned< Type > &)
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
static const pointConstraints & New(const pointMesh &mesh)
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
const Mesh & mesh() const
Return mesh.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Generic dimensioned Type class.
const double e
Elementary charge.
Definition: doubleFloat.H:78
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
#define forAll(list, i)
Definition: UList.H:421
label patchi
const dimensionSet & dimensions() const
Return dimensions.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const word & name() const
Return name.
Definition: IOobject.H:260
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
Generic GeometricField class.
#define forAllReverse(list, i)
Definition: UList.H:424
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Traits class for primitives.
Definition: pTraits.H:50
error FatalError
List< label > labelList
A List of labels.
Definition: labelList.H:56
const Time & time() const
Return time.
Definition: IOobject.C:251
fvOptions constrain(rhoEqn)
void smooth(const GeometricField< Type, pointPatchField, pointMesh > &fld, const scalarField &edgeWeight, GeometricField< Type, pointPatchField, pointMesh > &newFld) const
Fully explicit smoothing of fields (not positions)
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:245
A class for managing temporary objects.
Definition: PtrList.H:118
label nPoints() const
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
unsigned int get(const label) const
Get value at index I.
Definition: PackedListI.H:972