Loading...
Searching...
No Matches
cyclicACMIFvPatch.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) 2013-2016 OpenFOAM Foundation
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 "cyclicACMIFvPatch.H"
30#include "fvMesh.H"
31#include "transform.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
41}
42
43
44// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
45
47{
48 // Give AMI chance to update itself
49 bool updated = cyclicACMIPolyPatch_.updateAreas();
50
51 if (!cyclicACMIPolyPatch_.owner())
52 {
53 return updated;
54 }
55
56 if (updated || !cyclicACMIPolyPatch_.upToDate(areaTime_))
57 {
58 if (debug)
59 {
60 Pout<< "cyclicACMIFvPatch::updateAreas() : updating fv areas for "
61 << name() << " and " << this->nonOverlapPatch().name()
62 << endl;
63 }
64
65 const fvPatch& nonOverlapPatch = this->nonOverlapPatch();
66 const cyclicACMIFvPatch& nbrACMI = neighbPatch();
67 const fvPatch& nbrNonOverlapPatch = nbrACMI.nonOverlapPatch();
68
69 resetPatchAreas(*this);
71 resetPatchAreas(nbrACMI);
72 resetPatchAreas(nbrNonOverlapPatch);
73
74 updated = true;
75
76 // Mark my data to be up to date with ACMI polyPatch level
77 cyclicACMIPolyPatch_.setUpToDate(areaTime_);
78 }
79 return updated;
80}
81
82
83void Foam::cyclicACMIFvPatch::resetPatchAreas(const fvPatch& fvp) const
84{
85 const_cast<vectorField&>(fvp.Sf()) = fvp.patch().faceAreas();
86 const_cast<vectorField&>(fvp.Cf()) = fvp.patch().faceCentres();
87 const_cast<scalarField&>(fvp.magSf()) = mag(fvp.patch().faceAreas());
88
90 << fvp.patch().name() << " area:" << sum(fvp.magSf()) << endl;
91}
92
93
95{
96 if (coupled())
97 {
98 const cyclicACMIFvPatch& nbrPatch = neighbFvPatch();
99 const scalarField deltas(nf() & coupledFvPatch::delta());
100
101 // These deltas are of the cyclic part alone - they are
102 // not affected by the amount of overlap with the nonOverlapPatch
103 scalarField nbrDeltas
104 (
106 (
107 nbrPatch.nf() & nbrPatch.coupledFvPatch::delta()
108 )
109 );
110
111 const scalar tol = cyclicACMIPolyPatch::tolerance();
112
113
114 forAll(deltas, facei)
115 {
116 scalar di = mag(deltas[facei]);
117 scalar dni = mag(nbrDeltas[facei]);
118
119 if (dni < tol)
120 {
121 // Avoid zero weights on disconnected faces. This value
122 // will be weighted with the (zero) face area so will not
123 // influence calculations.
124 w[facei] = 1.0;
125 }
126 else
127 {
128 w[facei] = dni/(di + dni);
129 }
130 }
131 }
132 else
133 {
134 // Behave as uncoupled patch
136 }
137}
138
139
140// * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
141
143(
144 const polyPatch& patch,
145 const fvBoundaryMesh& bm
146)
147:
148 coupledFvPatch(patch, bm),
149 cyclicACMILduInterface(),
150 cyclicACMIPolyPatch_(refCast<const cyclicACMIPolyPatch>(patch)),
151 areaTime_
152 (
153 IOobject
154 (
155 "areaTime",
156 boundaryMesh().mesh().pointsInstance(),
157 boundaryMesh().mesh(),
158 IOobject::NO_READ,
159 IOobject::NO_WRITE,
160 IOobject::NO_REGISTER
161 ),
162 dimensionedScalar("time", dimTime, -GREAT)
163 )
164{
165 areaTime_.eventNo() = -1;
166}
167
168
169// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
170
171// void Foam::cyclicACMIFvPatch::newInternalProcFaces
172// (
173// label& newFaces,
174// label& newProcFaces
175// ) const
176// {
177// const List<labelList>& addSourceFaces =
178// cyclicACMIPolyPatch_.AMI().srcAddress();
179//
180// const scalarField& fMask = cyclicACMIPolyPatch_.srcMask();
181//
182// // Add new faces as many weights for AMI
183// forAll (addSourceFaces, faceI)
184// {
185// if (fMask[faceI] > cyclicACMIPolyPatch_.tolerance_)
186// {
187// const labelList& nbrFaceIs = addSourceFaces[faceI];
188//
189// forAll (nbrFaceIs, j)
190// {
191// label nbrFaceI = nbrFaceIs[j];
192//
193// if (nbrFaceI < neighbPatch().size())
194// {
195// // local faces
196// newFaces++;
197// }
198// else
199// {
200// // Proc faces
201// newProcFaces++;
202// }
203// }
204// }
205// }
206// }
207
208
209// Foam::refPtr<Foam::labelListList>
210// Foam::cyclicACMIFvPatch::mapCollocatedFaces() const
211// {
212// const scalarField& fMask = cyclicACMIPolyPatch_.srcMask();
213// const labelListList& srcFaces = cyclicACMIPolyPatch_.AMI().srcAddress();
214// labelListList dOverFaces;
215//
216// dOverFaces.setSize(srcFaces.size());
217// forAll (dOverFaces, faceI)
218// {
219// if (fMask[faceI] > cyclicACMIPolyPatch_.tolerance_)
220// {
221// dOverFaces[faceI].setSize(srcFaces[faceI].size());
222//
223// forAll (dOverFaces[faceI], subFaceI)
224// {
225// dOverFaces[faceI][subFaceI] = srcFaces[faceI][subFaceI];
226// }
227// }
228// }
229// return refPtr<labelListList>(new labelListList(dOverFaces));
230// }
231
234{
235 return Pstream::parRun() || (this->size() && neighbFvPatch().size());
236}
237
238
240{
241 if (coupled())
242 {
243 const cyclicACMIFvPatch& nbrPatch = neighbFvPatch();
244
245 const vectorField patchD(coupledFvPatch::delta());
246
247 vectorField nbrPatchD(interpolate(nbrPatch.coupledFvPatch::delta()));
248
249 auto tpdv = tmp<vectorField>::New(patchD.size());
250 vectorField& pdv = tpdv.ref();
251
252 // do the transformation if necessary
253 if (parallel())
254 {
255 forAll(patchD, facei)
256 {
257 const vector& ddi = patchD[facei];
258 const vector& dni = nbrPatchD[facei];
259
260 pdv[facei] = ddi - dni;
261 }
262 }
263 else
264 {
265 forAll(patchD, facei)
266 {
267 const vector& ddi = patchD[facei];
268 const vector& dni = nbrPatchD[facei];
269
270 pdv[facei] = ddi - transform(forwardT()[0], dni);
271 }
272 }
273
274 return tpdv;
275 }
276 else
277 {
278 return coupledFvPatch::delta();
279 }
280}
281
282
284(
285 const labelUList& internalData
286) const
287{
288 return patchInternalField(internalData);
289}
290
291
293(
294 const labelUList& internalData,
295 const labelUList& faceCells
296) const
298 auto tpfld = tmp<labelField>::New(this->size());
299 patchInternalField(internalData, faceCells, tpfld.ref());
300 return tpfld;
301}
302
303
305(
306 const Pstream::commsTypes commsType,
307 const labelUList& iF
308) const
309{
310 return neighbFvPatch().patchInternalField(iF);
311}
312
313
315{
316 // Update local and higher level areas
317 const bool updated = updateAreas();
318
319 // If anything changed update the mesh flux
320 if (cyclicACMIPolyPatch_.owner() && updated)
321 {
322 if (debug)
323 {
324 Pout<< "cyclicACMIFvPatch::movePoints() : areas updated for "
325 << name() << "; updating mesh flux now"
326 << endl;
327 }
328
329 // Scale the mesh flux
330
331 const fvPatch& nonOverlapPatch = this->nonOverlapPatch();
332 const cyclicACMIFvPatch& nbrACMI = neighbPatch();
333 const fvPatch& nbrNonOverlapPatch = nbrACMI.nonOverlapPatch();
334
335 const labelListList& newSrcAddr = AMI().srcAddress();
336 const labelListList& newTgtAddr = AMI().tgtAddress();
337
338 const fvMesh& mesh = boundaryMesh().mesh();
339 surfaceScalarField& meshPhi = const_cast<fvMesh&>(mesh).setPhi().ref();
340 surfaceScalarField::Boundary& meshPhiBf = meshPhi.boundaryFieldRef();
341
342 // Note: phip and phiNonOverlap will be different sizes if new faces
343 // have been added
344 scalarField& phip = meshPhiBf[cyclicACMIPolyPatch_.index()];
345 scalarField& phiNonOverlapp =
346 meshPhiBf[nonOverlapPatch.patch().index()];
347
348 const auto& points = mesh.points();
349
350 forAll(phip, facei)
351 {
352 if (newSrcAddr[facei].empty())
353 {
354 // AMI patch with no connection to other coupled faces
355 phip[facei] = 0.0;
356 }
357 else
358 {
359 // Scale the mesh flux according to the area fraction
360 const face& fAMI = cyclicACMIPolyPatch_[facei];
361
362 // Note: using raw point locations to calculate the geometric
363 // area - faces areas are currently scaled (decoupled from
364 // mesh points)
365 const scalar geomArea = fAMI.mag(points);
366 phip[facei] *= magSf()[facei]/geomArea;
367 }
368 }
369
370 forAll(phiNonOverlapp, facei)
371 {
372 const scalar w = 1.0 - cyclicACMIPolyPatch_.srcMask()[facei];
373 phiNonOverlapp[facei] *= w;
374 }
375
376 const cyclicACMIPolyPatch& nbrPatch = nbrACMI.cyclicACMIPatch();
377 scalarField& nbrPhip = meshPhiBf[nbrPatch.index()];
378 scalarField& nbrPhiNonOverlapp =
379 meshPhiBf[nbrNonOverlapPatch.patch().index()];
380
381 forAll(nbrPhip, facei)
382 {
383 if (newTgtAddr[facei].empty())
384 {
385 nbrPhip[facei] = 0.0;
386 }
387 else
388 {
389 const face& fAMI = nbrPatch[facei];
390
391 // Note: using raw point locations to calculate the geometric
392 // area - faces areas are currently scaled (decoupled from
393 // mesh points)
394 const scalar geomArea = fAMI.mag(points);
395 nbrPhip[facei] *= nbrACMI.magSf()[facei]/geomArea;
396 }
397 }
398
399 forAll(nbrPhiNonOverlapp, facei)
400 {
401 const scalar w = 1.0 - cyclicACMIPolyPatch_.tgtMask()[facei];
402 nbrPhiNonOverlapp[facei] *= w;
403 }
404 }
405}
406
407
408// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
commsTypes
Communications types.
Definition UPstream.H:81
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
const bMesh & mesh() const
virtual const labelUList & faceCells() const
Return faceCell addressing.
coupledFvPatch(const polyPatch &patch, const fvBoundaryMesh &bm)
Construct from polyPatch.
virtual tmp< vectorField > delta() const =0
Return delta (P to N) vectors across coupled patch.
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI).
const cyclicACMIFvPatch & neighbFvPatch() const
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
void resetPatchAreas(const fvPatch &fvp) const
Helper function to reset the FV patch areas from the primitive patch.
void makeWeights(scalarField &) const
Make patch weighting factors.
virtual void movePoints()
Correct patches after moving points.
virtual bool coupled() const
Return true if this patch is coupled.
virtual bool parallel() const
Are the cyclic planes parallel.
const cyclicACMIPolyPatch & cyclicACMIPatch() const
Return local reference cast into the cyclic patch.
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
cyclicACMIFvPatch(const polyPatch &patch, const fvBoundaryMesh &bm)
Construct from polyPatch.
virtual const fvPatch & nonOverlapPatch() const
Return non-overlapping fvPatch.
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Return the values of the given internal data adjacent to the interface as a field.
virtual const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
virtual const cyclicACMIFvPatch & neighbPatch() const
Return neighbour fvPatch.
virtual bool updateAreas() const
Update the AMI and patch areas. Return true if anything updated.
tmp< Field< Type > > interpolate(const Field< Type > &fld) const
Interpolate (make sure to have uptodate areas).
virtual const tensorField & forwardT() const
Return face transformation tensor.
cyclicACMILduInterface() noexcept=default
Default construct.
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI).
bool upToDate(const regIOobject &) const
Return true if given object is up to date with *this.
static scalar tolerance()
Overlap tolerance.
virtual bool updateAreas() const
Update the AMI and patch areas. Return true if anything.
virtual bool owner() const
Does this side own the patch?
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
scalar mag(const UList< point > &p) const
Magnitude of face area.
Definition faceI.H:105
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual label size() const
Patch size is the number of faces, but can be overloaded.
Definition fvPatch.H:242
virtual const word & name() const
Return name.
Definition fvPatch.H:210
const polyPatch & patch() const noexcept
Return the polyPatch.
Definition fvPatch.H:202
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition fvPatch.C:157
const scalarField & magSf() const
Return face area magnitudes, like the fvMesh::magSf() method.
Definition fvPatch.C:131
tmp< vectorField > nf() const
Return face unit normals, like the fvMesh::unitSf() method Same as unitSf().
Definition fvPatch.C:143
const vectorField & Cf() const
Return face centres.
Definition fvPatch.C:113
void patchInternalField(const UList< Type > &internalData, const labelUList &addressing, UList< Type > &pfld) const
Extract internal field next to patch using specified addressing.
const vectorField & Sf() const
Return face area vectors, like the fvMesh::Sf() method.
Definition fvPatch.C:125
label index() const noexcept
The index of this patch in the boundaryMesh.
const word & name() const noexcept
The patch name.
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const vectorField::subField faceAreas() const
Return face normals.
Definition polyPatch.C:326
const vectorField::subField faceCentres() const
Return face centres.
Definition polyPatch.C:320
label eventNo() const noexcept
Event number at last update.
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
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
bool coupled
dynamicFvMesh & mesh
auto & name
const pointField & points
#define DebugPout
Report an information message using Foam::Pout.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
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)
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition curveTools.C:75
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.
3D tensor transformation operations.