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 
55 {
56  return
58  || !this->boundaryMesh().mesh().time().processorCase();
59 }
60 
61 
63 {
64  if (coupled())
65  {
66  const cyclicAMIFvPatch& nbrPatch = neighbFvPatch();
67 
68  const scalarField deltas(nf() & coupledFvPatch::delta());
69 
70  tmp<scalarField> tnbrDeltas;
71  if (applyLowWeightCorrection())
72  {
73  tnbrDeltas =
75  (
76  nbrPatch.nf() & nbrPatch.coupledFvPatch::delta(),
77  scalarField(this->size(), 1.0)
78  );
79  }
80  else
81  {
82  tnbrDeltas =
83  interpolate(nbrPatch.nf() & nbrPatch.coupledFvPatch::delta());
84  }
85 
86  const scalarField& nbrDeltas = tnbrDeltas();
87 
88  forAll(deltas, facei)
89  {
90  // Note use of mag
91  scalar di = mag(deltas[facei]);
92  scalar dni = mag(nbrDeltas[facei]);
93 
94  w[facei] = dni/(di + dni);
95  }
96  }
97  else
98  {
99  // Behave as uncoupled patch
101  }
102 }
103 
104 
106 {
107  // Apply correction to default coeffs
108 }
109 
110 
112 {
113  // Apply correction to default coeffs
114  //coeffs = Zero;
115 }
116 
117 
119 {
120  // Apply correction to default vectors
121  //vecs = Zero;
122 }
123 
124 
126 {
127  const cyclicAMIFvPatch& nbrPatch = neighbFvPatch();
128 
129  if (coupled())
130  {
131  const vectorField patchD(coupledFvPatch::delta());
132 
133  tmp<vectorField> tnbrPatchD;
134  if (applyLowWeightCorrection())
135  {
136  tnbrPatchD =
138  (
139  nbrPatch.coupledFvPatch::delta(),
140  vectorField(this->size(), Zero)
141  );
142  }
143  else
144  {
145  tnbrPatchD = interpolate(nbrPatch.coupledFvPatch::delta());
146  }
147 
148  const vectorField& nbrPatchD = tnbrPatchD();
149 
150  auto tpdv = tmp<vectorField>::New(patchD.size());
151  vectorField& pdv = tpdv.ref();
152 
153  // do the transformation if necessary
154  if (parallel())
155  {
156  forAll(patchD, facei)
157  {
158  const vector& ddi = patchD[facei];
159  const vector& dni = nbrPatchD[facei];
160 
161  pdv[facei] = ddi - dni;
162  }
163  }
164  else
165  {
166  forAll(patchD, facei)
167  {
168  const vector& ddi = patchD[facei];
169  const vector& dni = nbrPatchD[facei];
170 
171  pdv[facei] = ddi - transform(forwardT()[0], dni);
172  }
173  }
174 
175  return tpdv;
176  }
177  else
178  {
179  return coupledFvPatch::delta();
180  }
181 }
182 
183 
185 (
186  const labelUList& internalData
187 ) const
188 {
189  return patchInternalField(internalData);
190 }
191 
192 
194 (
195  const Pstream::commsTypes commsType,
196  const labelUList& iF
197 ) const
198 {
199  return neighbFvPatch().patchInternalField(iF);
200 }
201 
202 
204 {
205  if (!owner() || !cyclicAMIPolyPatch_.createAMIFaces())
206  {
207  // Only manipulating patch face areas and mesh motion flux if the AMI
208  // creates additional faces
209  return;
210  }
211 
212  // Update face data based on values set by the AMI manipulations
213  const_cast<vectorField&>(Sf()) = cyclicAMIPolyPatch_.faceAreas();
214  const_cast<vectorField&>(Cf()) = cyclicAMIPolyPatch_.faceCentres();
215  const_cast<scalarField&>(magSf()) = mag(Sf());
216 
217  const cyclicAMIFvPatch& nbr = neighbPatch();
218  const_cast<vectorField&>(nbr.Sf()) = nbr.cyclicAMIPatch().faceAreas();
219  const_cast<vectorField&>(nbr.Cf()) = nbr.cyclicAMIPatch().faceCentres();
220  const_cast<scalarField&>(nbr.magSf()) = mag(nbr.Sf());
221 
222 
223  // Set consitent mesh motion flux
224  // TODO: currently maps src mesh flux to tgt - update to
225  // src = src + mapped(tgt) and tgt = tgt + mapped(src)?
226 
227  const fvMesh& mesh = boundaryMesh().mesh();
228  surfaceScalarField& meshPhi = const_cast<fvMesh&>(mesh).setPhi();
229  surfaceScalarField::Boundary& meshPhiBf = meshPhi.boundaryFieldRef();
230 
231  if (cyclicAMIPolyPatch_.owner())
232  {
233  scalarField& phip = meshPhiBf[patch().index()];
234  forAll(phip, facei)
235  {
236  const face& f = cyclicAMIPolyPatch_.localFaces()[facei];
237 
238  // Note: using raw point locations to calculate the geometric
239  // area - faces areas are currently scaled by the AMI weights
240  // (decoupled from mesh points)
241  const scalar geomArea = f.mag(cyclicAMIPolyPatch_.localPoints());
242 
243  const scalar scaledArea = magSf()[facei];
244  phip[facei] *= scaledArea/geomArea;
245  }
246 
247  scalarField srcMeshPhi(phip);
248  if (Pstream::parRun())
249  {
250  AMI().srcMap().distribute(srcMeshPhi);
251  }
252 
253  const labelListList& tgtToSrcAddr = AMI().tgtAddress();
254  scalarField& nbrPhip = meshPhiBf[nbr.index()];
255 
256  forAll(tgtToSrcAddr, tgtFacei)
257  {
258  // Note: now have 1-to-1 mapping so tgtToSrcAddr[tgtFacei] is size 1
259  const label srcFacei = tgtToSrcAddr[tgtFacei][0];
260  nbrPhip[tgtFacei] = -srcMeshPhi[srcFacei];
261  }
262 
263  DebugInfo
264  << "patch:" << patch().name()
265  << " sum(area):" << gSum(magSf())
266  << " min(mag(faceAreas):" << gMin(magSf())
267  << " sum(meshPhi):" << gSum(phip) << nl
268  << " sum(nbrMeshPhi):" << gSum(nbrPhip) << nl
269  << endl;
270  }
271 }
272 
273 
274 // ************************************************************************* //
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:194
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:125
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:111
Foam::UPstream::parRun
static bool & parRun()
Test if this a parallel run, or allow modify access.
Definition: UPstream.H:434
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:105
Foam::fvBoundaryMesh::mesh
const fvMesh & mesh() const
Return the mesh reference.
Definition: fvBoundaryMesh.H: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:118
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: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:54
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:83
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Time.H
Foam::cyclicAMIFvPatch::interfaceInternalField
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Definition: cyclicAMIFvPatch.C:185
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:69
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:359
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:62
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::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
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::TimePaths::processorCase
bool processorCase() const
Return true if this is a processor case.
Definition: TimePathsI.H:36
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:203