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-------------------------------------------------------------------------------
11License
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"
30#include "coupledPolyPatch.H"
32#include "transform.H"
33#include "faMesh.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
41}
42
43const Foam::scalar Foam::cyclicFaPatch::matchTol_ = 1e-3;
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
49(
50 const word& name,
51 const dictionary& dict,
52 const label index,
53 const faBoundaryMesh& bm,
54 const word& patchType
55)
56:
57 coupledFaPatch(name, dict, index, bm, patchType)
58{}
59
60
61
62// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63
64void Foam::cyclicFaPatch::calcTransforms()
65{
66 const label sizeby2 = this->size()/2;
67 if (size() > 0)
68 {
69 pointField half0Ctrs(sizeby2);
70 pointField half1Ctrs(sizeby2);
71 for (label edgei=0; edgei < sizeby2; ++edgei)
72 {
73 half0Ctrs[edgei] = this->edgeCentres()[edgei];
74 half1Ctrs[edgei] = this->edgeCentres()[edgei + sizeby2];
75 }
76
77 vectorField half0Normals(sizeby2);
78 vectorField half1Normals(sizeby2);
79
80 const vectorField eN(edgeNormals()*magEdgeLengths());
81
82 scalar maxMatchError = 0;
83 label errorEdge = -1;
84
85 for (label edgei = 0; edgei < sizeby2; ++edgei)
86 {
87 half0Normals[edgei] = eN[edgei];
88 half1Normals[edgei] = eN[edgei + sizeby2];
89
90 scalar magLe = mag(half0Normals[edgei]);
91 scalar nbrMagLe = mag(half1Normals[edgei]);
92 scalar avLe = 0.5*(magLe + nbrMagLe);
93
94 if (magLe < ROOTVSMALL && nbrMagLe < ROOTVSMALL)
95 {
96 // Undetermined normal. Use dummy normal to force separation
97 // check. (note use of sqrt(VSMALL) since that is how mag
98 // scales)
99 half0Normals[edgei] = point(1, 0, 0);
100 half1Normals[edgei] = half0Normals[edgei];
101 }
102 else if (mag(magLe - nbrMagLe)/avLe > matchTol_)
103 {
104 // Error in area matching. Find largest error
105 maxMatchError =
106 Foam::max(maxMatchError, mag(magLe - nbrMagLe)/avLe);
107 errorEdge = edgei;
108 }
109 else
110 {
111 half0Normals[edgei] /= magLe;
112 half1Normals[edgei] /= nbrMagLe;
113 }
114 }
115
116 // Check for error in edge matching
117 if (maxMatchError > matchTol_)
118 {
119 label nbrEdgei = errorEdge + sizeby2;
120 scalar magLe = mag(half0Normals[errorEdge]);
121 scalar nbrMagLe = mag(half1Normals[errorEdge]);
122 scalar avLe = 0.5*(magLe + nbrMagLe);
123
125 << "edge " << errorEdge
126 << " area does not match neighbour "
127 << nbrEdgei << " by "
128 << 100*mag(magLe - nbrMagLe)/avLe
129 << "% -- possible edge ordering problem." << endl
130 << "patch:" << name()
131 << " my area:" << magLe
132 << " neighbour area:" << nbrMagLe
133 << " matching tolerance:" << matchTol_
134 << endl
135 << "Mesh edge:" << start() + errorEdge
136 << endl
137 << "Neighbour edge:" << start() + nbrEdgei
138 << endl
139 << "Other errors also exist, only the largest is reported. "
140 << "Please rerun with cyclic debug flag set"
141 << " for more information." << exit(FatalError);
142 }
143
144 // Calculate transformation tensors
145 calcTransformTensors
146 (
147 half0Ctrs,
148 half1Ctrs,
149 half0Normals,
150 half1Normals
151 );
152
153 // Check transformation tensors
154 if (!parallel())
155 {
156 if (forwardT().size() > 1 || reverseT().size() > 1)
157 {
159 << "Transformation tensor is not constant for the cyclic "
160 << "patch. Please reconsider your setup and definition of "
161 << "cyclic boundaries." << endl;
162 }
163 }
164 }
165}
166
167
169{
170 const scalarField& magL = magEdgeLengths();
171
172 const scalarField deltas(edgeNormals() & faPatch::delta());
173 const label sizeby2 = deltas.size()/2;
174
175 scalar maxMatchError = 0;
176 label errorEdge = -1;
177
178 for (label edgei = 0; edgei < sizeby2; ++edgei)
179 {
180 scalar avL = 0.5*(magL[edgei] + magL[edgei + sizeby2]);
181
182 if
183 (
184 mag(magL[edgei] - magL[edgei + sizeby2])/avL
185 > matchTol_
186 )
187 {
188 // Found error. Look for largest matching error
189 maxMatchError =
191 (
192 maxMatchError,
193 mag(magL[edgei] - magL[edgei + sizeby2])/avL
194 );
195
196 errorEdge = edgei;
197 }
198
199 scalar di = deltas[edgei];
200 scalar dni = deltas[edgei + sizeby2];
201
202 w[edgei] = dni/(di + dni);
203 w[edgei + sizeby2] = 1 - w[edgei];
204 }
205
206 // Check for error in matching
207 if (maxMatchError > matchTol_)
208 {
209 scalar avL = 0.5*(magL[errorEdge] + magL[errorEdge + sizeby2]);
210
212 << "edge " << errorEdge << " and " << errorEdge + sizeby2
213 << " areas do not match by "
214 << 100*mag(magL[errorEdge] - magL[errorEdge + sizeby2])/avL
215 << "% -- possible edge ordering problem." << nl
216 << "Cyclic area match tolerance = "
217 << matchTol_ << " patch: " << name()
218 << abort(FatalError);
219 }
220}
221
222
224{
225 const scalarField deltas(edgeNormals() & faPatch::delta());
226 const label sizeby2 = deltas.size()/2;
227
228 for (label edgei = 0; edgei < sizeby2; ++edgei)
229 {
230 scalar di = deltas[edgei];
231 scalar dni = deltas[edgei + sizeby2];
232
233 dc[edgei] = 1.0/(di + dni);
234 dc[edgei + sizeby2] = dc[edgei];
235 }
236}
237
238
240{
242}
243
244
246{
248 calcTransforms();
249}
250
251
253(
254 PstreamBuffers& pBufs,
255 const pointField& p
256)
257{
259}
260
261
263(
264 PstreamBuffers& pBufs,
265 const pointField& p
266)
267{
268 faPatch::movePoints(pBufs, p);
269 calcTransforms();
270}
271
272
274{
275 const vectorField patchD(faPatch::delta());
276 const label sizeby2 = patchD.size()/2;
277
278 auto tpdv = tmp<vectorField>::New(patchD.size());
279 auto& pdv = tpdv.ref();
280
281 // Do the transformation if necessary
282 if (parallel())
283 {
284 for (label edgei = 0; edgei < sizeby2; ++edgei)
285 {
286 const vector& ddi = patchD[edgei];
287 vector dni = patchD[edgei + sizeby2];
288
289 pdv[edgei] = ddi - dni;
290 pdv[edgei + sizeby2] = -pdv[edgei];
291 }
292 }
293 else
294 {
295 for (label edgei = 0; edgei < sizeby2; ++edgei)
296 {
297 const vector& ddi = patchD[edgei];
298 vector dni = patchD[edgei + sizeby2];
299
300 pdv[edgei] = ddi - transform(forwardT()[0], dni);
301 pdv[edgei + sizeby2] = -transform(reverseT()[0], pdv[edgei]);
302 }
303 }
304
305 return tpdv;
306}
307
308
310(
311 const labelUList& internalData
312) const
313{
314 return patchInternalField(internalData);
315}
316
317
319(
320 const labelUList& internalData,
321 const labelUList& edgeFaces
322) const
323{
324 return patchInternalField(internalData, edgeFaces);
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// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
commsTypes
Types of communications.
Definition: UPstream.H:67
coupledFaPatch is an abstract base class for patches that couple regions of the computational domain ...
Cyclic-plane patch.
Definition: cyclicFaPatch.H:61
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
void makeWeights(scalarField &) const
Make patch weighting factors.
void makeDeltaCoeffs(scalarField &) const
Make patch face - neighbour cell distances.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
static const scalar matchTol_
Relative tolerance (for geometric matching). Is factor of.
Definition: cyclicFaPatch.H:72
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Finite area boundary mesh.
void calcGeometry()
Calculate the geometry for the patches.
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition: faPatch.H:78
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: faPatch.H:123
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: faPatch.C:494
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
virtual void initMovePoints()
Initialise the patches for moving points.
Definition: fvPatch.C:188
void movePoints()
Update for new mesh geometry.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
vector point
Point is a vector.
Definition: point.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Field< vector > vectorField
Specialisation of Field<T> for vector.
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
volScalarField & e
Definition: createFields.H:11
3D tensor transformation operations.