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 "transform.H"
33 #include "surfaceFields.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(cyclicAMIFvPatch, 0);
40  addToRunTimeSelectionTable(fvPatch, cyclicAMIFvPatch, polyPatch);
42  (
43  fvPatch,
44  cyclicAMIFvPatch,
45  polyPatch,
46  cyclicPeriodicAMI
47  );
48 }
49 
50 
51 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
52 
54 {
55  return Pstream::parRun() || (this->size() && neighbFvPatch().size());
56 }
57 
58 
60 {
61  if (coupled())
62  {
63  const cyclicAMIFvPatch& nbrPatch = neighbFvPatch();
64 
65  const scalarField deltas(nf() & coupledFvPatch::delta());
66 
67  tmp<scalarField> tnbrDeltas;
68  if (applyLowWeightCorrection())
69  {
70  tnbrDeltas =
72  (
73  nbrPatch.nf() & nbrPatch.coupledFvPatch::delta(),
74  scalarField(this->size(), 1.0)
75  );
76  }
77  else
78  {
79  tnbrDeltas =
80  interpolate(nbrPatch.nf() & nbrPatch.coupledFvPatch::delta());
81  }
82 
83  const scalarField& nbrDeltas = tnbrDeltas();
84 
85  forAll(deltas, facei)
86  {
87  // Note use of mag
88  scalar di = mag(deltas[facei]);
89  scalar dni = mag(nbrDeltas[facei]);
90 
91  w[facei] = dni/(di + dni);
92  }
93  }
94  else
95  {
96  // Behave as uncoupled patch
98  }
99 }
100 
101 
103 {
104  // Apply correction to default coeffs
105 }
106 
107 
109 {
110  // Apply correction to default coeffs
111  //coeffs = Zero;
112 }
113 
114 
116 {
117  // Apply correction to default vectors
118  //vecs = Zero;
119 }
120 
121 
123 {
124  const cyclicAMIFvPatch& nbrPatch = neighbFvPatch();
125 
126  if (coupled())
127  {
128  const vectorField patchD(coupledFvPatch::delta());
129 
130  tmp<vectorField> tnbrPatchD;
131  if (applyLowWeightCorrection())
132  {
133  tnbrPatchD =
135  (
136  nbrPatch.coupledFvPatch::delta(),
137  vectorField(this->size(), Zero)
138  );
139  }
140  else
141  {
142  tnbrPatchD = interpolate(nbrPatch.coupledFvPatch::delta());
143  }
144 
145  const vectorField& nbrPatchD = tnbrPatchD();
146 
147  auto tpdv = tmp<vectorField>::New(patchD.size());
148  vectorField& pdv = tpdv.ref();
149 
150  // do the transformation if necessary
151  if (parallel())
152  {
153  forAll(patchD, facei)
154  {
155  const vector& ddi = patchD[facei];
156  const vector& dni = nbrPatchD[facei];
157 
158  pdv[facei] = ddi - dni;
159  }
160  }
161  else
162  {
163  forAll(patchD, facei)
164  {
165  const vector& ddi = patchD[facei];
166  const vector& dni = nbrPatchD[facei];
167 
168  pdv[facei] = ddi - transform(forwardT()[0], dni);
169  }
170  }
171 
172  return tpdv;
173  }
174  else
175  {
176  return coupledFvPatch::delta();
177  }
178 }
179 
180 
182 (
183  const labelUList& internalData
184 ) const
185 {
186  return patchInternalField(internalData);
187 }
188 
189 
191 (
192  const Pstream::commsTypes commsType,
193  const labelUList& iF
194 ) const
195 {
196  return neighbFvPatch().patchInternalField(iF);
197 }
198 
199 
201 {
202  if (!owner() || !cyclicAMIPolyPatch_.createAMIFaces())
203  {
204  // Only manipulating patch face areas and mesh motion flux if the AMI
205  // creates additional faces
206  return;
207  }
208 
209  // Update face data based on values set by the AMI manipulations
210  const_cast<vectorField&>(Sf()) = cyclicAMIPolyPatch_.faceAreas();
211  const_cast<vectorField&>(Cf()) = cyclicAMIPolyPatch_.faceCentres();
212  const_cast<scalarField&>(magSf()) = mag(Sf());
213 
214  const cyclicAMIFvPatch& nbr = neighbPatch();
215  const_cast<vectorField&>(nbr.Sf()) = nbr.cyclicAMIPatch().faceAreas();
216  const_cast<vectorField&>(nbr.Cf()) = nbr.cyclicAMIPatch().faceCentres();
217  const_cast<scalarField&>(nbr.magSf()) = mag(nbr.Sf());
218 
219 
220  // Set consitent mesh motion flux
221  // TODO: currently maps src mesh flux to tgt - update to
222  // src = src + mapped(tgt) and tgt = tgt + mapped(src)?
223 
224  const fvMesh& mesh = boundaryMesh().mesh();
225  surfaceScalarField& meshPhi = const_cast<fvMesh&>(mesh).setPhi();
226  surfaceScalarField::Boundary& meshPhiBf = meshPhi.boundaryFieldRef();
227 
228  if (cyclicAMIPolyPatch_.owner())
229  {
230  scalarField& phip = meshPhiBf[patch().index()];
231  forAll(phip, facei)
232  {
233  const face& f = cyclicAMIPolyPatch_.localFaces()[facei];
234 
235  // Note: using raw point locations to calculate the geometric
236  // area - faces areas are currently scaled by the AMI weights
237  // (decoupled from mesh points)
238  const scalar geomArea = f.mag(cyclicAMIPolyPatch_.localPoints());
239 
240  const scalar scaledArea = magSf()[facei];
241  phip[facei] *= scaledArea/geomArea;
242  }
243 
244  scalarField srcMeshPhi(phip);
245  if (Pstream::parRun())
246  {
247  AMI().srcMap().distribute(srcMeshPhi);
248  }
249 
250  const labelListList& tgtToSrcAddr = AMI().tgtAddress();
251  scalarField& nbrPhip = meshPhiBf[nbr.index()];
252 
253  forAll(tgtToSrcAddr, tgtFacei)
254  {
255  // Note: now have 1-to-1 mapping so tgtToSrcAddr[tgtFacei] is size 1
256  const label srcFacei = tgtToSrcAddr[tgtFacei][0];
257  nbrPhip[tgtFacei] = -srcMeshPhi[srcFacei];
258  }
259 
260  DebugInfo
261  << "patch:" << patch().name()
262  << " sum(area):" << gSum(magSf())
263  << " min(mag(faceAreas):" << gMin(magSf())
264  << " sum(meshPhi):" << gSum(phip) << nl
265  << " sum(nbrMeshPhi):" << gSum(nbrPhip) << nl
266  << endl;
267  }
268 }
269 
270 
271 // ************************************************************************* //
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:144
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:191
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
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:122
Foam::coupledFvPatch::delta
virtual tmp< vectorField > delta() const =0
Return delta (P to N) vectors across coupled patch.
Definition: coupledFvPatch.C:46
Foam::cyclicAMIFvPatch::makeNonOrthoDeltaCoeffs
virtual void makeNonOrthoDeltaCoeffs(scalarField &) const
Correct patch non-ortho deltaCoeffs.
Definition: cyclicAMIFvPatch.C:108
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:415
Foam::cyclicAMIFvPatch::neighbFvPatch
const cyclicAMIFvPatch & neighbFvPatch() const
Definition: cyclicAMIFvPatch.H:161
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:350
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:102
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
Foam::cyclicAMIFvPatch::makeNonOrthoCorrVectors
virtual void makeNonOrthoCorrVectors(vectorField &) const
Correct patch non-ortho correction vectors.
Definition: cyclicAMIFvPatch.C:115
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::size
virtual label size() const
Return size.
Definition: fvPatch.H:175
Foam::fvPatch::nf
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:138
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:113
Foam::cyclicAMIFvPatch::coupled
virtual bool coupled() const
Definition: cyclicAMIFvPatch.C:53
Foam::fvPatch::index
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:193
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:84
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::cyclicAMIFvPatch::interfaceInternalField
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Definition: cyclicAMIFvPatch.C:182
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:66
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:354
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::polyPatch::faceCentres
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:313
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::List< labelList >
Foam::cyclicAMIFvPatch::makeWeights
void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: cyclicAMIFvPatch.C:59
Foam::polyPatch::faceAreas
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:319
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::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:105
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:53
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:150
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::fvPatch::makeWeights
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: fvPatch.C:164
cyclicAMIFvPatch.H
Foam::cyclicAMIFvPatch::movePoints
virtual void movePoints()
Correct patches after moving points.
Definition: cyclicAMIFvPatch.C:200