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-2020 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 
136 
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  return patchInternalField(internalData, faceCells);
232 }
233 
234 
236 (
237  const Pstream::commsTypes commsType,
238  const labelUList& iF
239 ) const
240 {
241  return neighbFvPatch().patchInternalField(iF);
242 }
243 
244 
246 {
247  if (!owner() || !cyclicAMIPolyPatch_.createAMIFaces())
248  {
249  // Only manipulating patch face areas and mesh motion flux if the AMI
250  // creates additional faces
251  return;
252  }
253 
254  // Update face data based on values set by the AMI manipulations
255  const_cast<vectorField&>(Sf()) = cyclicAMIPolyPatch_.faceAreas();
256  const_cast<vectorField&>(Cf()) = cyclicAMIPolyPatch_.faceCentres();
257  const_cast<scalarField&>(magSf()) = mag(Sf());
258 
259  const cyclicAMIFvPatch& nbr = neighbPatch();
260  const_cast<vectorField&>(nbr.Sf()) = nbr.cyclicAMIPatch().faceAreas();
261  const_cast<vectorField&>(nbr.Cf()) = nbr.cyclicAMIPatch().faceCentres();
262  const_cast<scalarField&>(nbr.magSf()) = mag(nbr.Sf());
263 
264 
265  // Set consitent mesh motion flux
266  // TODO: currently maps src mesh flux to tgt - update to
267  // src = src + mapped(tgt) and tgt = tgt + mapped(src)?
268 
269  const fvMesh& mesh = boundaryMesh().mesh();
270  surfaceScalarField& meshPhi = const_cast<fvMesh&>(mesh).setPhi();
271  surfaceScalarField::Boundary& meshPhiBf = meshPhi.boundaryFieldRef();
272 
273  if (cyclicAMIPolyPatch_.owner())
274  {
275  scalarField& phip = meshPhiBf[patch().index()];
276  forAll(phip, facei)
277  {
278  const face& f = cyclicAMIPolyPatch_.localFaces()[facei];
279 
280  // Note: using raw point locations to calculate the geometric
281  // area - faces areas are currently scaled by the AMI weights
282  // (decoupled from mesh points)
283  const scalar geomArea = f.mag(cyclicAMIPolyPatch_.localPoints());
284 
285  const scalar scaledArea = magSf()[facei];
286  phip[facei] *= scaledArea/geomArea;
287  }
288 
289  scalarField srcMeshPhi(phip);
290  if (AMI().distributed())
291  {
292  AMI().srcMap().distribute(srcMeshPhi);
293  }
294 
295  const labelListList& tgtToSrcAddr = AMI().tgtAddress();
296  scalarField& nbrPhip = meshPhiBf[nbr.index()];
297 
298  forAll(tgtToSrcAddr, tgtFacei)
299  {
300  // Note: now have 1-to-1 mapping so tgtToSrcAddr[tgtFacei] is size 1
301  const label srcFacei = tgtToSrcAddr[tgtFacei][0];
302  nbrPhip[tgtFacei] = -srcMeshPhi[srcFacei];
303  }
304 
305  DebugInfo
306  << "patch:" << patch().name()
307  << " sum(area):" << gSum(magSf())
308  << " min(mag(faceAreas):" << gMin(magSf())
309  << " sum(meshPhi):" << gSum(phip) << nl
310  << " sum(nbrMeshPhi):" << gSum(nbrPhip) << nl
311  << endl;
312  }
313 }
314 
315 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
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::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::cyclicAMIFvPatch::internalFieldTransfer
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
Definition: cyclicAMIFvPatch.C:236
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::cyclicAMIFvPatch::delta
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
Definition: cyclicAMIFvPatch.C:157
Foam::coupledFvPatch::delta
virtual tmp< vectorField > delta() const =0
Return delta (P to N) vectors across coupled patch.
Definition: coupledFvPatch.C:46
Foam::TimePaths::processorCase
bool processorCase() const noexcept
Return true if this is a processor case.
Definition: TimePathsI.H:36
Foam::cyclicAMIFvPatch::makeNonOrthoDeltaCoeffs
virtual void makeNonOrthoDeltaCoeffs(scalarField &) const
Correct patch non-ortho deltaCoeffs.
Definition: cyclicAMIFvPatch.C:143
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::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::cyclicAMIFvPatch::makeDeltaCoeffs
virtual void makeDeltaCoeffs(scalarField &) const
Correct patch deltaCoeffs.
Definition: cyclicAMIFvPatch.C:137
Foam::fvBoundaryMesh::mesh
const fvMesh & mesh() const noexcept
Return the mesh reference.
Definition: fvBoundaryMesh.H:103
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Foam::cyclicAMIFvPatch::makeNonOrthoCorrVectors
virtual void makeNonOrthoCorrVectors(vectorField &) const
Correct patch non-ortho correction vectors.
Definition: cyclicAMIFvPatch.C:150
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::fvPatch::boundaryMesh
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:203
Foam::fvPatch::nf
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:144
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::Field< scalar >
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
Foam::fvPatch::Cf
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:119
Foam::cyclicAMIFvPatch::coupled
virtual bool coupled() const
Definition: cyclicAMIFvPatch.C:86
Foam::fvPatch::index
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:197
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::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:69
Time.H
Foam::cyclicAMIFvPatch::interfaceInternalField
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Definition: cyclicAMIFvPatch.C:217
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::polyPatch::faceCentres
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:321
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< labelList >
Foam::cyclicAMIFvPatch::makeWeights
void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: cyclicAMIFvPatch.C:94
Foam::polyPatch::faceAreas
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:327
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::UList< label >
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
Foam::cyclicAMIFvPatch::cyclicAMIPatch
const cyclicAMIPolyPatch & cyclicAMIPatch() const
Return local reference cast into the cyclic patch.
Definition: cyclicAMIFvPatch.H:141
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::cyclicAMIFvPatch
Cyclic patch for Arbitrary Mesh Interface (AMI)
Definition: cyclicAMIFvPatch.H:54
coupled
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
transform.H
3D tensor transformation operations.
Foam::fvPatch::magSf
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:156
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
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
cyclicAMIFvPatch.H
Foam::cyclicAMIFvPatch::movePoints
virtual void movePoints()
Correct patches after moving points.
Definition: cyclicAMIFvPatch.C:245