cyclicAMIFvPatch.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2023 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 "cyclicAMIFvPatch.H"
31 #include "fvMesh.H"
32 #include "Time.H"
33 #include "transform.H"
34 #include "surfaceFields.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(cyclicAMIFvPatch, 0);
41  addToRunTimeSelectionTable(fvPatch, cyclicAMIFvPatch, polyPatch);
43  (
44  fvPatch,
45  cyclicAMIFvPatch,
46  polyPatch,
47  cyclicPeriodicAMI
48  );
49 }
50 
51 
52 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
53 
54 // void Foam::cyclicAMIFvPatch::newInternalProcFaces
55 // (
56 // label& newFaces,
57 // label& newProcFaces
58 // ) const
59 // {
60 // const labelListList& addSourceFaces = AMI().srcAddress();
61 //
62 // // Add new faces as many weights for AMI
63 // forAll (addSourceFaces, faceI)
64 // {
65 // const labelList& nbrFaceIs = addSourceFaces[faceI];
66 //
67 // forAll (nbrFaceIs, j)
68 // {
69 // label nbrFaceI = nbrFaceIs[j];
70 //
71 // if (nbrFaceI < neighbPatch().size())
72 // {
73 // // local faces
74 // newFaces++;
75 // }
76 // else
77 // {
78 // // Proc faces
79 // newProcFaces++;
80 // }
81 // }
82 // }
83 // }
84 
85 
87 {
88  return
90  || !this->boundaryMesh().mesh().time().processorCase();
91 }
92 
93 
95 {
96  if (coupled())
97  {
98  const cyclicAMIFvPatch& nbrPatch = neighbFvPatch();
99 
100  const scalarField deltas(nf() & coupledFvPatch::delta());
101 
102  tmp<scalarField> tnbrDeltas;
103  if (applyLowWeightCorrection())
104  {
105  tnbrDeltas =
107  (
108  nbrPatch.nf() & nbrPatch.coupledFvPatch::delta(),
109  scalarField(this->size(), 1.0)
110  );
111  }
112  else
113  {
114  tnbrDeltas =
115  interpolate(nbrPatch.nf() & nbrPatch.coupledFvPatch::delta());
116  }
117 
118  const scalarField& nbrDeltas = tnbrDeltas();
119 
120  forAll(deltas, facei)
121  {
122  // Note use of mag
123  scalar di = mag(deltas[facei]);
124  scalar dni = mag(nbrDeltas[facei]);
125 
126  w[facei] = dni/(di + dni);
127  }
128  }
129  else
130  {
131  // Behave as uncoupled patch
133  }
134 }
135 
138 {
139  // Apply correction to default coeffs
140 }
141 
142 
144 {
145  // Apply correction to default coeffs
146  //coeffs = Zero;
147 }
148 
149 
151 {
152  // Apply correction to default vectors
153  //vecs = Zero;
154 }
155 
156 
158 {
159  const cyclicAMIFvPatch& nbrPatch = neighbFvPatch();
160 
161  if (coupled())
162  {
163  const vectorField patchD(coupledFvPatch::delta());
164 
165  tmp<vectorField> tnbrPatchD;
166  if (applyLowWeightCorrection())
167  {
168  tnbrPatchD =
170  (
171  nbrPatch.coupledFvPatch::delta(),
172  vectorField(this->size(), Zero)
173  );
174  }
175  else
176  {
177  tnbrPatchD = interpolate(nbrPatch.coupledFvPatch::delta());
178  }
179 
180  const vectorField& nbrPatchD = tnbrPatchD();
181 
182  auto tpdv = tmp<vectorField>::New(patchD.size());
183  vectorField& pdv = tpdv.ref();
184 
185  // do the transformation if necessary
186  if (parallel())
187  {
188  forAll(patchD, facei)
189  {
190  const vector& ddi = patchD[facei];
191  const vector& dni = nbrPatchD[facei];
192 
193  pdv[facei] = ddi - dni;
194  }
195  }
196  else
197  {
198  forAll(patchD, facei)
199  {
200  const vector& ddi = patchD[facei];
201  const vector& dni = nbrPatchD[facei];
202 
203  pdv[facei] = ddi - transform(forwardT()[0], dni);
204  }
205  }
206 
207  return tpdv;
208  }
209  else
210  {
211  return coupledFvPatch::delta();
212  }
213 }
214 
215 
217 (
218  const labelUList& internalData
219 ) const
220 {
221  return patchInternalField(internalData);
222 }
223 
224 
226 (
227  const labelUList& internalData,
228  const labelUList& faceCells
229 ) const
230 {
231  auto tpfld = tmp<labelField>::New();
232  patchInternalField(internalData, faceCells, tpfld.ref());
233  return tpfld;
234 }
235 
236 
238 (
239  const Pstream::commsTypes commsType,
240  const labelUList& iF
241 ) const
242 {
243  return neighbFvPatch().patchInternalField(iF);
244 }
245 
246 
248 {
249  if (!owner() || !cyclicAMIPolyPatch_.createAMIFaces())
250  {
251  // Only manipulating patch face areas and mesh motion flux if the AMI
252  // creates additional faces
253  return;
254  }
255 
256  // Update face data based on values set by the AMI manipulations
257  const_cast<vectorField&>(Sf()) = cyclicAMIPolyPatch_.faceAreas();
258  const_cast<vectorField&>(Cf()) = cyclicAMIPolyPatch_.faceCentres();
259  const_cast<scalarField&>(magSf()) = mag(Sf());
260 
261  const cyclicAMIFvPatch& nbr = neighbPatch();
262  const_cast<vectorField&>(nbr.Sf()) = nbr.cyclicAMIPatch().faceAreas();
263  const_cast<vectorField&>(nbr.Cf()) = nbr.cyclicAMIPatch().faceCentres();
264  const_cast<scalarField&>(nbr.magSf()) = mag(nbr.Sf());
265 
266 
267  // Set consitent mesh motion flux
268  // TODO: currently maps src mesh flux to tgt - update to
269  // src = src + mapped(tgt) and tgt = tgt + mapped(src)?
270 
271  const fvMesh& mesh = boundaryMesh().mesh();
272  surfaceScalarField& meshPhi = const_cast<fvMesh&>(mesh).setPhi().ref();
273  surfaceScalarField::Boundary& meshPhiBf = meshPhi.boundaryFieldRef();
274 
275  if (cyclicAMIPolyPatch_.owner())
276  {
277  scalarField& phip = meshPhiBf[patch().index()];
278  forAll(phip, facei)
279  {
280  const face& f = cyclicAMIPolyPatch_.localFaces()[facei];
281 
282  // Note: using raw point locations to calculate the geometric
283  // area - faces areas are currently scaled by the AMI weights
284  // (decoupled from mesh points)
285  const scalar geomArea = f.mag(cyclicAMIPolyPatch_.localPoints());
286 
287  const scalar scaledArea = magSf()[facei];
288  phip[facei] *= scaledArea/geomArea;
289  }
290 
291  scalarField srcMeshPhi(phip);
292  if (AMI().distributed())
293  {
294  AMI().srcMap().distribute(srcMeshPhi);
295  }
296 
297  const labelListList& tgtToSrcAddr = AMI().tgtAddress();
298  scalarField& nbrPhip = meshPhiBf[nbr.index()];
299 
300  forAll(tgtToSrcAddr, tgtFacei)
301  {
302  // Note: now have 1-to-1 mapping so tgtToSrcAddr[tgtFacei] is size 1
303  const label srcFacei = tgtToSrcAddr[tgtFacei][0];
304  nbrPhip[tgtFacei] = -srcMeshPhi[srcFacei];
305  }
306 
307  DebugInfo
308  << "patch:" << patch().name()
309  << " sum(area):" << gSum(magSf())
310  << " min(mag(faceAreas):" << gMin(magSf())
311  << " sum(meshPhi):" << gSum(phip) << nl
312  << " sum(nbrMeshPhi):" << gSum(nbrPhip) << nl
313  << endl;
314  }
315 }
316 
317 // ************************************************************************* //
Foam::surfaceFields.
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
commsTypes
Communications types.
Definition: UPstream.H:72
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Return the values of the given internal data adjacent to the interface as a field.
Type gMin(const FieldField< Field, Type > &f)
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
tmp< vectorField > nf() const
Return face unit normals, like the fvMesh::unitSf() method Same as unitSf().
Definition: fvPatch.C:143
Cyclic patch for Arbitrary Mesh Interface (AMI)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
const fvMesh & mesh() const noexcept
Return the mesh reference.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
virtual void makeDeltaCoeffs(scalarField &) const
Correct patch deltaCoeffs.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
bool processorCase() const noexcept
True if this is a processor case.
Definition: TimePathsI.H:52
virtual void makeNonOrthoCorrVectors(vectorField &) const
Correct patch non-ortho correction vectors.
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
3D tensor transformation operations.
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: fvPatch.C:157
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
Vector< scalar > vector
Definition: vector.H:57
#define DebugInfo
Report an information message using Foam::Info.
const fvBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition: fvPatch.H:260
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
virtual void makeNonOrthoDeltaCoeffs(scalarField &) const
Correct patch non-ortho deltaCoeffs.
virtual bool coupled() const
Return true if this patch is coupled. This is equivalent to the coupledPolyPatch::coupled() if parall...
virtual tmp< vectorField > delta() const =0
Return delta (P to N) vectors across coupled patch.
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:29
const std::string patch
OpenFOAM patch number as a std::string.
virtual void movePoints()
Correct patches after moving points.
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
void makeWeights(scalarField &) const
Make patch weighting factors.
Field< vector > vectorField
Specialisation of Field<T> for vector.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
bool coupled
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127