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 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 "cyclicACMIFvPatch.H"
30 #include "fvMesh.H"
31 #include "transform.H"
32 #include "surfaceFields.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(cyclicACMIFvPatch, 0);
40  addToRunTimeSelectionTable(fvPatch, cyclicACMIFvPatch, polyPatch);
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 
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 
89  DebugPout
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),
150  cyclicACMIPolyPatch_(refCast<const cyclicACMIPolyPatch>(patch)),
151  areaTime_
152  (
153  IOobject
154  (
155  "areaTime",
156  boundaryMesh().mesh().pointsInstance(),
157  boundaryMesh().mesh(),
160  false
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 
232 
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
297 {
298  return patchInternalField(internalData, faceCells);
299 }
300 
301 
303 (
304  const Pstream::commsTypes commsType,
305  const labelUList& iF
306 ) const
307 {
308  return neighbFvPatch().patchInternalField(iF);
309 }
310 
311 
313 {
314  if (!cyclicACMIPolyPatch_.owner())
315  {
316  return;
317  }
318 
319 
320  if (!cyclicACMIPolyPatch_.upToDate(areaTime_))
321  {
322  if (debug)
323  {
324  Pout<< "cyclicACMIFvPatch::movePoints() : updating fv areas for "
325  << name() << " and " << this->nonOverlapPatch().name()
326  << endl;
327  }
328 
329 
330  // Set the patch face areas to be consistent with the changes made
331  // at the polyPatch level
332 
333  const fvPatch& nonOverlapPatch = this->nonOverlapPatch();
334  const cyclicACMIFvPatch& nbrACMI = neighbPatch();
335  const fvPatch& nbrNonOverlapPatch = nbrACMI.nonOverlapPatch();
336 
337  resetPatchAreas(*this);
338  resetPatchAreas(nonOverlapPatch);
339  resetPatchAreas(nbrACMI);
340  resetPatchAreas(nbrNonOverlapPatch);
341 
342  // Scale the mesh flux
343 
344  const labelListList& newSrcAddr = AMI().srcAddress();
345  const labelListList& newTgtAddr = AMI().tgtAddress();
346 
347  const fvMesh& mesh = boundaryMesh().mesh();
348  surfaceScalarField& meshPhi = const_cast<fvMesh&>(mesh).setPhi();
349  surfaceScalarField::Boundary& meshPhiBf = meshPhi.boundaryFieldRef();
350 
351  // Note: phip and phiNonOverlap will be different sizes if new faces
352  // have been added
353  scalarField& phip = meshPhiBf[cyclicACMIPolyPatch_.index()];
354  scalarField& phiNonOverlapp =
355  meshPhiBf[nonOverlapPatch.patch().index()];
356 
357  const auto& points = mesh.points();
358 
359  forAll(phip, facei)
360  {
361  if (newSrcAddr[facei].empty())
362  {
363  // AMI patch with no connection to other coupled faces
364  phip[facei] = 0.0;
365  }
366  else
367  {
368  // Scale the mesh flux according to the area fraction
369  const face& fAMI = cyclicACMIPolyPatch_[facei];
370 
371  // Note: using raw point locations to calculate the geometric
372  // area - faces areas are currently scaled (decoupled from
373  // mesh points)
374  const scalar geomArea = fAMI.mag(points);
375  phip[facei] *= magSf()[facei]/geomArea;
376  }
377  }
378 
379  forAll(phiNonOverlapp, facei)
380  {
381  const scalar w = 1.0 - cyclicACMIPolyPatch_.srcMask()[facei];
382  phiNonOverlapp[facei] *= w;
383  }
384 
385  const cyclicACMIPolyPatch& nbrPatch = nbrACMI.cyclicACMIPatch();
386  scalarField& nbrPhip = meshPhiBf[nbrPatch.index()];
387  scalarField& nbrPhiNonOverlapp =
388  meshPhiBf[nbrNonOverlapPatch.patch().index()];
389 
390  forAll(nbrPhip, facei)
391  {
392  if (newTgtAddr[facei].empty())
393  {
394  nbrPhip[facei] = 0.0;
395  }
396  else
397  {
398  const face& fAMI = nbrPatch[facei];
399 
400  // Note: using raw point locations to calculate the geometric
401  // area - faces areas are currently scaled (decoupled from
402  // mesh points)
403  const scalar geomArea = fAMI.mag(points);
404  nbrPhip[facei] *= nbrACMI.magSf()[facei]/geomArea;
405  }
406  }
407 
408  forAll(nbrPhiNonOverlapp, facei)
409  {
410  const scalar w = 1.0 - cyclicACMIPolyPatch_.tgtMask()[facei];
411  nbrPhiNonOverlapp[facei] *= w;
412  }
413 
414  // Mark my data to be up to date with ACMI polyPatch level
415  cyclicACMIPolyPatch_.setUpToDate(areaTime_);
416  }
417 }
418 
419 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::boundaryMesh::mesh
const bMesh & mesh() const
Definition: boundaryMesh.H:206
Foam::fvPatch::Sf
const vectorField & Sf() const
Return face area vectors.
Definition: fvPatch.C:150
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::cyclicAMIPolyPatch::owner
virtual bool owner() const
Does this side own the patch?
Definition: cyclicAMIPolyPatch.C:964
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::face::mag
scalar mag(const UList< point > &p) const
Magnitude of face area.
Definition: faceI.H:112
Foam::cyclicACMIPolyPatch::upToDate
bool upToDate(const regIOobject &) const
Return true if given object is up to date with *this.
Definition: cyclicACMIPolyPatch.C:142
Foam::coupledFvPatch::delta
virtual tmp< vectorField > delta() const =0
Return delta (P to N) vectors across coupled patch.
Definition: coupledFvPatch.C:46
Foam::cyclicACMILduInterface
An abstract base class for cyclic ACMI coupled interfaces.
Definition: cyclicACMILduInterface.H:51
Foam::cyclicACMIPolyPatch::tolerance
static scalar tolerance()
Overlap tolerance.
Definition: cyclicACMIPolyPatchI.H:65
Foam::cyclicACMIPolyPatch::setUpToDate
void setUpToDate(regIOobject &) const
Set object up to date with *this.
Definition: cyclicACMIPolyPatch.C:151
Foam::fvc::meshPhi
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:36
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
surfaceFields.H
Foam::surfaceFields.
Foam::fvPatch::name
virtual const word & name() const
Return name.
Definition: fvPatch.H:167
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Foam::cyclicACMIFvPatch::delta
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
Definition: cyclicACMIFvPatch.C:239
Foam::cyclicACMIFvPatch::coupled
virtual bool coupled() const
Return true if this patch is coupled.
Definition: cyclicACMIFvPatch.C:233
Foam::cyclicACMIFvPatch::internalFieldTransfer
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
Definition: cyclicACMIFvPatch.C:303
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:57
Foam::fvPatch::nf
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:144
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::Field< vector >
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::cyclicACMIFvPatch::neighbPatch
virtual const cyclicACMIFvPatch & neighbPatch() const
Return neighbour fvPatch.
Definition: cyclicACMIFvPatch.H:128
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::fvPatch::Cf
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:119
Foam::cyclicACMIFvPatch::movePoints
virtual void movePoints()
Correct patches after moving points.
Definition: cyclicACMIFvPatch.C:312
Foam::cyclicACMIFvPatch::nonOverlapPatch
virtual const fvPatch & nonOverlapPatch() const
Return non-overlapping fvPatch.
Definition: cyclicACMIFvPatch.H:143
Foam::cyclicACMIFvPatch::updateAreas
virtual bool updateAreas() const
Update the AMI and patch areas. Return true if anything updated.
Definition: cyclicACMIFvPatch.C:46
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::patchIdentifier::index
label index() const noexcept
The index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:147
DebugPout
#define DebugPout
Report an information message using Foam::Pout.
Definition: messageStream.H:393
Foam::cyclicACMIFvPatch::cyclicACMIFvPatch
cyclicACMIFvPatch(const polyPatch &patch, const fvBoundaryMesh &bm)
Construct from polyPatch.
Definition: cyclicACMIFvPatch.C:143
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:69
Foam::polyPatch::faceCentres
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:321
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::cyclicACMIFvPatch::resetPatchAreas
void resetPatchAreas(const fvPatch &fvp) const
Helper function to reset the FV patch areas from the primitive patch.
Definition: cyclicACMIFvPatch.C:83
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< labelList >
Foam::polyPatch::faceAreas
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:327
Foam::cyclicACMIPolyPatch::updateAreas
virtual bool updateAreas() const
Update the AMI and patch areas. Return true if anything.
Definition: cyclicACMIPolyPatch.C:46
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::cyclicACMIFvPatch::cyclicACMIPatch
const cyclicACMIPolyPatch & cyclicACMIPatch() const
Return local reference cast into the cyclic patch.
Definition: cyclicACMIFvPatch.H:111
Foam::UList< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::cyclicACMIFvPatch::interfaceInternalField
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Definition: cyclicACMIFvPatch.C:284
Foam::fvPatch::patch
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:161
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
cyclicACMIFvPatch.H
Foam::patchIdentifier::name
const word & name() const noexcept
The patch name.
Definition: patchIdentifier.H:135
Foam::coupledFvPatch
An abstract base class for patches that couple regions of the computational domain e....
Definition: coupledFvPatch.H:53
coupled
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
Foam::cyclicACMIPolyPatch
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI).
Definition: cyclicACMIPolyPatch.H:78
transform.H
3D tensor transformation operations.
Foam::fvPatch::magSf
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:156
Foam::cyclicACMIFvPatch::makeWeights
void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: cyclicACMIFvPatch.C:94
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
Foam::fvPatch::makeWeights
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: fvPatch.C:170
Foam::cyclicACMIFvPatch
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI)
Definition: cyclicACMIFvPatch.H:54