cyclicFaPatch.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) 2016-2017 Wikki Ltd
9  Copyright (C) 2019-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "cyclicFaPatch.H"
31 #include "transform.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(cyclicFaPatch, 0);
38  addToRunTimeSelectionTable(faPatch, cyclicFaPatch, dictionary);
39 }
40 
41 const Foam::scalar Foam::cyclicFaPatch::matchTol_ = 1e-3;
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const word& name,
49  const dictionary& dict,
50  const label index,
51  const faBoundaryMesh& bm,
52  const word& patchType
53 )
54 :
55  coupledFaPatch(name, dict, index, bm, patchType)
56 {}
57 
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
62 void Foam::cyclicFaPatch::calcTransforms()
63 {
64  const label sizeby2 = this->size()/2;
65  if (size() > 0)
66  {
67  pointField half0Ctrs(sizeby2);
68  pointField half1Ctrs(sizeby2);
69  for (label edgei=0; edgei < sizeby2; ++edgei)
70  {
71  half0Ctrs[edgei] = this->edgeCentres()[edgei];
72  half1Ctrs[edgei] = this->edgeCentres()[edgei + sizeby2];
73  }
74 
75  vectorField half0Normals(sizeby2);
76  vectorField half1Normals(sizeby2);
77 
78  const vectorField eN(edgeNormals()*magEdgeLengths());
79 
80  scalar maxMatchError = 0;
81  label errorEdge = -1;
82 
83  for (label edgei = 0; edgei < sizeby2; ++edgei)
84  {
85  half0Normals[edgei] = eN[edgei];
86  half1Normals[edgei] = eN[edgei + sizeby2];
87 
88  scalar magLe = mag(half0Normals[edgei]);
89  scalar nbrMagLe = mag(half1Normals[edgei]);
90  scalar avLe = 0.5*(magLe + nbrMagLe);
91 
92  if (magLe < ROOTVSMALL && nbrMagLe < ROOTVSMALL)
93  {
94  // Undetermined normal. Use dummy normal to force separation
95  // check. (note use of sqrt(VSMALL) since that is how mag
96  // scales)
97  half0Normals[edgei] = point(1, 0, 0);
98  half1Normals[edgei] = half0Normals[edgei];
99  }
100  else if (mag(magLe - nbrMagLe)/avLe > matchTol_)
101  {
102  // Error in area matching. Find largest error
103  maxMatchError =
104  Foam::max(maxMatchError, mag(magLe - nbrMagLe)/avLe);
105  errorEdge = edgei;
106  }
107  else
108  {
109  half0Normals[edgei] /= magLe;
110  half1Normals[edgei] /= nbrMagLe;
111  }
112  }
113 
114  // Check for error in edge matching
115  if (maxMatchError > matchTol_)
116  {
117  label nbrEdgei = errorEdge + sizeby2;
118  scalar magLe = mag(half0Normals[errorEdge]);
119  scalar nbrMagLe = mag(half1Normals[errorEdge]);
120  scalar avLe = 0.5*(magLe + nbrMagLe);
121 
123  << "edge " << errorEdge
124  << " area does not match neighbour "
125  << nbrEdgei << " by "
126  << 100*mag(magLe - nbrMagLe)/avLe
127  << "% -- possible edge ordering problem." << endl
128  << "patch:" << name()
129  << " my area:" << magLe
130  << " neighbour area:" << nbrMagLe
131  << " matching tolerance:" << matchTol_
132  << endl
133  << "Mesh edge:" << start() + errorEdge
134  << endl
135  << "Neighbour edge:" << start() + nbrEdgei
136  << endl
137  << "Other errors also exist, only the largest is reported. "
138  << "Please rerun with cyclic debug flag set"
139  << " for more information." << exit(FatalError);
140  }
141 
142  // Calculate transformation tensors
143  calcTransformTensors
144  (
145  half0Ctrs,
146  half1Ctrs,
147  half0Normals,
148  half1Normals
149  );
150 
151  // Check transformation tensors
152  if (!parallel())
153  {
154  if (forwardT().size() > 1 || reverseT().size() > 1)
155  {
157  << "Transformation tensor is not constant for the cyclic "
158  << "patch. Please reconsider your setup and definition of "
159  << "cyclic boundaries." << endl;
160  }
161  }
162  }
163 }
164 
165 
167 {
168  const scalarField& magL = magEdgeLengths();
169 
170  const scalarField deltas(edgeNormals() & coupledFaPatch::delta());
171  const label sizeby2 = deltas.size()/2;
172 
173  scalar maxMatchError = 0;
174  label errorEdge = -1;
175 
176  for (label edgei = 0; edgei < sizeby2; ++edgei)
177  {
178  scalar avL = 0.5*(magL[edgei] + magL[edgei + sizeby2]);
179 
180  if
181  (
182  mag(magL[edgei] - magL[edgei + sizeby2])/avL
183  > matchTol_
184  )
185  {
186  // Found error. Look for largest matching error
187  maxMatchError =
188  Foam::max
189  (
190  maxMatchError,
191  mag(magL[edgei] - magL[edgei + sizeby2])/avL
192  );
193 
194  errorEdge = edgei;
195  }
196 
197  scalar di = deltas[edgei];
198  scalar dni = deltas[edgei + sizeby2];
199 
200  w[edgei] = dni/(di + dni);
201  w[edgei + sizeby2] = 1 - w[edgei];
202  }
203 
204  // Check for error in matching
205  if (maxMatchError > matchTol_)
206  {
207  scalar avL = 0.5*(magL[errorEdge] + magL[errorEdge + sizeby2]);
208 
210  << "edge " << errorEdge << " and " << errorEdge + sizeby2
211  << " areas do not match by "
212  << 100*mag(magL[errorEdge] - magL[errorEdge + sizeby2])/avL
213  << "% -- possible edge ordering problem." << nl
214  << "Cyclic area match tolerance = "
215  << matchTol_ << " patch: " << name()
216  << abort(FatalError);
217  }
218 }
219 
220 
222 {
223  const scalarField deltas(edgeNormals() & coupledFaPatch::delta());
224  const label sizeby2 = deltas.size()/2;
225 
226  for (label edgei = 0; edgei < sizeby2; ++edgei)
227  {
228  scalar di = deltas[edgei];
229  scalar dni = deltas[edgei + sizeby2];
231  dc[edgei] = 1.0/(di + dni);
232  dc[edgei + sizeby2] = dc[edgei];
233  }
234 }
235 
238 {
239  faPatch::initGeometry(pBufs);
240 }
241 
242 
244 {
245  faPatch::calcGeometry(pBufs);
246  calcTransforms();
247 }
248 
249 
251 (
252  PstreamBuffers& pBufs,
253  const pointField& p
254 )
255 {
256  faPatch::initMovePoints(pBufs, p);
257 }
258 
259 
261 (
262  PstreamBuffers& pBufs,
263  const pointField& p
264 )
265 {
266  faPatch::movePoints(pBufs, p);
267  calcTransforms();
268 }
269 
270 
272 {
273  const vectorField patchD(coupledFaPatch::delta());
274  const label sizeby2 = patchD.size()/2;
275 
276  auto tpdv = tmp<vectorField>::New(patchD.size());
277  auto& pdv = tpdv.ref();
278 
279  // Do the transformation if necessary
280  if (parallel())
281  {
282  for (label edgei = 0; edgei < sizeby2; ++edgei)
283  {
284  const vector& ddi = patchD[edgei];
285  const vector& dni = patchD[edgei + sizeby2];
286 
287  pdv[edgei] = ddi - dni;
288  pdv[edgei + sizeby2] = -pdv[edgei];
289  }
290  }
291  else
292  {
293  for (label edgei = 0; edgei < sizeby2; ++edgei)
294  {
295  const vector& ddi = patchD[edgei];
296  const vector& dni = patchD[edgei + sizeby2];
297 
298  pdv[edgei] = ddi - transform(forwardT()[0], dni);
299  pdv[edgei + sizeby2] = -transform(reverseT()[0], pdv[edgei]);
300  }
301  }
302 
303  return tpdv;
304 }
305 
306 
308 (
309  const labelUList& internalData
310 ) const
311 {
312  return patchInternalField(internalData);
313 }
314 
315 
317 (
318  const labelUList& internalData,
319  const labelUList& edgeFaces
320 ) const
321 {
322  auto tpfld = tmp<labelField>::New();
323  patchInternalField(internalData, edgeFaces, tpfld.ref());
324  return tpfld;
325 }
326 
327 
329 (
330  const Pstream::commsTypes commsType,
331  const labelUList& interfaceData
332 ) const
333 {
334  auto tpnf = tmp<labelField>::New(this->size());
335  auto& pnf = tpnf.ref();
336 
337  const label sizeby2 = this->size()/2;
338 
339  for (label edgei=0; edgei < sizeby2; ++edgei)
340  {
341  pnf[edgei] = interfaceData[edgei + sizeby2];
342  pnf[edgei + sizeby2] = interfaceData[edgei];
343  }
344 
345  return tpnf;
346 }
347 
348 
350 (
351  const Pstream::commsTypes commsType,
352  const labelUList& iF
353 ) const
354 {
355  return internalFieldTransfer(commsType, iF, this->faceCells());
356 }
357 
358 
360 (
361  const Pstream::commsTypes commsType,
362  const labelUList& iF,
363  const labelUList& edgeCells
364 ) const
365 {
366  auto tpnf = tmp<labelField>::New(this->size());
367  auto& pnf = tpnf.ref();
368 
369  const label sizeby2 = this->size()/2;
370 
371  for (label edgei=0; edgei < sizeby2; ++edgei)
372  {
373  pnf[edgei] = iF[edgeCells[edgei + sizeby2]];
374  pnf[edgei + sizeby2] = iF[edgeCells[edgei]];
375  }
376 
377  return tpnf;
378 }
379 
380 
381 // ************************************************************************* //
dictionary dict
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Return the values of the given internal data adjacent to the interface as a field.
static const scalar matchTol_
Relative tolerance (for geometric matching). Is factor of.
Definition: cyclicFaPatch.H:71
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
commsTypes
Communications types.
Definition: UPstream.H:77
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patches after moving points.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
void makeDeltaCoeffs(scalarField &) const
Make patch face - neighbour cell distances.
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
Definition: faPatch.H:157
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Definition: faPatch.H:151
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
virtual tmp< vectorField > delta() const =0
Return delta (P to N) vectors across coupled patch.
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ const Field< vector > edgeCentres(faMeshTools::flattenEdgeField(aMesh.edgeCentres(), true))
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
3D tensor transformation operations.
void makeWeights(scalarField &) const
Make patch weighting factors.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
Vector< scalar > vector
Definition: vector.H:57
errorManip< error > abort(error &err)
Definition: errorManip.H:139
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
defineTypeNameAndDebug(combustionModel, 0)
Buffers for inter-processor communications streams (UOPstream, UIPstream).
cyclicFaPatch(const word &name, const dictionary &dict, const label index, const faBoundaryMesh &bm, const word &patchType)
Construct from dictionary.
Definition: cyclicFaPatch.C:40
vector point
Point is a vector.
Definition: point.H:37
coupledFaPatch is an abstract base class for patches that couple regions of the computational domain ...
Finite area boundary mesh.
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual tmp< labelField > transfer(const Pstream::commsTypes commsType, const labelUList &interfaceData) const
Transfer and return neighbour field.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patch after moving points.
Definition: faPatch.C:528
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: faPatch.H:145
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)