Loading...
Searching...
No Matches
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-2025 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"
31#include "transform.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
40
41const 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
62void 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 =
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{
224 lPN = scalar(1)/lPN;
225}
226
227
229{
230 const scalarField deltas(edgeNormals() & coupledFaPatch::delta());
231 const label sizeby2 = deltas.size()/2;
232
233 for (label edgei = 0; edgei < sizeby2; ++edgei)
234 {
235 scalar di = deltas[edgei];
236 scalar dni = deltas[edgei + sizeby2];
238 dc[edgei] = 1.0/(di + dni);
239 dc[edgei + sizeby2] = dc[edgei];
240 }
241}
242
248
249
251{
253 calcTransforms();
254}
255
256
258(
259 PstreamBuffers& pBufs,
261)
262{
264}
265
266
268(
269 PstreamBuffers& pBufs,
270 const pointField& p
272{
273 faPatch::movePoints(pBufs, p);
274 calcTransforms();
275}
276
277
279{
280 const vectorField patchD(coupledFaPatch::delta());
281 const label sizeby2 = patchD.size()/2;
282
283 auto tpdv = tmp<vectorField>::New(patchD.size());
284 auto& pdv = tpdv.ref();
285
286 // Do the transformation if necessary
287 if (parallel())
288 {
289 for (label edgei = 0; edgei < sizeby2; ++edgei)
290 {
291 const vector& ddi = patchD[edgei];
292 const vector& dni = patchD[edgei + sizeby2];
293
294 pdv[edgei] = ddi - dni;
295 pdv[edgei + sizeby2] = -pdv[edgei];
296 }
297 }
298 else
299 {
300 for (label edgei = 0; edgei < sizeby2; ++edgei)
301 {
302 const vector& ddi = patchD[edgei];
303 const vector& dni = patchD[edgei + sizeby2];
304
305 pdv[edgei] = ddi - transform(forwardT()[0], dni);
306 pdv[edgei + sizeby2] = -transform(reverseT()[0], pdv[edgei]);
308 }
309
310 return tpdv;
311}
312
313
314Foam::tmp<Foam::labelField> Foam::cyclicFaPatch::interfaceInternalField
315(
316 const labelUList& internalData
317) const
318{
319 return patchInternalField(internalData);
320}
321
322
324(
325 const labelUList& internalData,
326 const labelUList& edgeFaces
327) const
329 auto tpfld = tmp<labelField>::New(this->size());
330 patchInternalField(internalData, edgeFaces, tpfld.ref());
331 return tpfld;
332}
333
334
336(
337 const Pstream::commsTypes commsType,
338 const labelUList& interfaceData
339) const
340{
341 auto tpnf = tmp<labelField>::New(this->size());
342 auto& pnf = tpnf.ref();
343
344 const label sizeby2 = this->size()/2;
345
346 for (label edgei=0; edgei < sizeby2; ++edgei)
347 {
348 pnf[edgei] = interfaceData[edgei + sizeby2];
349 pnf[edgei + sizeby2] = interfaceData[edgei];
350 }
351
352 return tpnf;
353}
354
355
357(
358 const Pstream::commsTypes commsType,
359 const labelUList& iF
360) const
361{
362 return internalFieldTransfer(commsType, iF, this->faceCells());
363}
364
365
367(
368 const Pstream::commsTypes commsType,
369 const labelUList& iF,
370 const labelUList& edgeCells
371) const
372{
373 auto tpnf = tmp<labelField>::New(this->size());
374 auto& pnf = tpnf.ref();
375
376 const label sizeby2 = this->size()/2;
377
378 for (label edgei=0; edgei < sizeby2; ++edgei)
379 {
380 pnf[edgei] = iF[edgeCells[edgei + sizeby2]];
381 pnf[edgei + sizeby2] = iF[edgeCells[edgei]];
382 }
383
384 return tpnf;
385}
386
387
388// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Info<< " "<< writer.output().name()<< nl;}{ const Field< vector > edgeCentres(faMeshTools::flattenEdgeField(aMesh.edgeCentres(), true))
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
commsTypes
Communications types.
Definition UPstream.H:81
coupledFaPatch is an abstract base class for patches that couple regions of the computational domain ...
bool parallel() const
Are the cyclic planes parallel.
virtual const labelUList & faceCells() const
Return faceCell addressing: lduInterface virtual function.
coupledFaPatch(const word &name, const labelUList &edgeLabels, const label index, const faBoundaryMesh &bm, const label nbrPolyPatchIndex, const word &patchType)
Construct from components.
virtual tmp< vectorField > delta() const =0
Return delta (P to N) vectors across coupled patch.
Cyclic-plane patch.
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
void makeWeights(scalarField &) const
Make patch weighting factors.
cyclicFaPatch(const word &name, const dictionary &dict, const label index, const faBoundaryMesh &bm, const word &patchType)
Construct from dictionary.
void makeLPN(scalarField &) const
Make patch geodesic distance between P and N.
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
void makeDeltaCoeffs(scalarField &) const
Make patch face - neighbour cell distances.
virtual tmp< labelField > transfer(const Pstream::commsTypes commsType, const labelUList &interfaceData) const
Transfer and return neighbour field.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
static const scalar matchTol_
Relative tolerance (for geometric matching). Is factor of.
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Return the values of the given internal data adjacent to the interface as a field.
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patches after moving points.
virtual const tensorField & forwardT() const
Return face transformation tensor.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Finite area boundary mesh, which is a faPatch list with registered IO, a reference to the associated ...
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition faPatch.H:76
virtual label size() const
Patch size is the number of edge labels, but can be overloaded.
Definition faPatch.H:392
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition faPatch.H:145
const scalarField & magEdgeLengths() const
Return edge length magnitudes, like the faMesh::magLe() method.
Definition faPatch.C:450
const scalarField & lPN() const
Return patch geodesic distance between P and N.
Definition faPatch.C:502
const labelUList & edgeFaces() const
Return edge-face addressing.
Definition faPatch.C:424
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patch after moving points.
Definition faPatch.C:547
tmp< vectorField > edgeNormals() const
Return edge unit normals, like the faMesh::unitLe() method.
Definition faPatch.C:456
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Definition faPatch.H:151
void patchInternalField(const UList< Type > &internalData, const labelUList &addressing, UList< Type > &pfld) const
Extract internal field next to patch using specified addressing.
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
Definition faPatch.H:157
label index() const noexcept
The index of this patch in the boundaryMesh.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
#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:40
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
volScalarField & e
3D tensor transformation operations.