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  const_cast<vectorField&>(fvp.Sf()) = fvp.patch().faceAreas();
49  const_cast<vectorField&>(fvp.Cf()) = fvp.patch().faceCentres();
50  const_cast<scalarField&>(fvp.magSf()) = mag(fvp.patch().faceAreas());
51 
52  DebugPout
53  << fvp.patch().name() << " area:" << sum(fvp.magSf()) << endl;
54 }
55 
56 
58 {
59  if (coupled())
60  {
61  const cyclicACMIFvPatch& nbrPatch = neighbFvPatch();
62  const scalarField deltas(nf() & coupledFvPatch::delta());
63 
64  // These deltas are of the cyclic part alone - they are
65  // not affected by the amount of overlap with the nonOverlapPatch
66  scalarField nbrDeltas
67  (
69  (
70  nbrPatch.nf() & nbrPatch.coupledFvPatch::delta()
71  )
72  );
73 
74  const scalar tol = cyclicACMIPolyPatch::tolerance();
75 
76 
77  forAll(deltas, facei)
78  {
79  scalar di = mag(deltas[facei]);
80  scalar dni = mag(nbrDeltas[facei]);
81 
82  if (dni < tol)
83  {
84  // Avoid zero weights on disconnected faces. This value
85  // will be weighted with the (zero) face area so will not
86  // influence calculations.
87  w[facei] = 1.0;
88  }
89  else
90  {
91  w[facei] = dni/(di + dni);
92  }
93  }
94  }
95  else
96  {
97  // Behave as uncoupled patch
99  }
100 }
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
106 {
107  return Pstream::parRun() || (this->size() && neighbFvPatch().size());
108 }
109 
110 
112 {
113  if (coupled())
114  {
115  const cyclicACMIFvPatch& nbrPatch = neighbFvPatch();
116 
117  const vectorField patchD(coupledFvPatch::delta());
118 
119  vectorField nbrPatchD(interpolate(nbrPatch.coupledFvPatch::delta()));
120 
121  auto tpdv = tmp<vectorField>::New(patchD.size());
122  vectorField& pdv = tpdv.ref();
123 
124  // do the transformation if necessary
125  if (parallel())
126  {
127  forAll(patchD, facei)
128  {
129  const vector& ddi = patchD[facei];
130  const vector& dni = nbrPatchD[facei];
131 
132  pdv[facei] = ddi - dni;
133  }
134  }
135  else
136  {
137  forAll(patchD, facei)
138  {
139  const vector& ddi = patchD[facei];
140  const vector& dni = nbrPatchD[facei];
141 
142  pdv[facei] = ddi - transform(forwardT()[0], dni);
143  }
144  }
145 
146  return tpdv;
147  }
148  else
149  {
150  return coupledFvPatch::delta();
151  }
152 }
153 
154 
156 (
157  const labelUList& internalData
158 ) const
159 {
160  return patchInternalField(internalData);
161 }
162 
163 
165 (
166  const Pstream::commsTypes commsType,
167  const labelUList& iF
168 ) const
169 {
170  return neighbFvPatch().patchInternalField(iF);
171 }
172 
173 
175 {
176  if (!cyclicACMIPolyPatch_.owner())
177  {
178  return;
179  }
180 
181  // Set the patch face areas to be consistent with the changes made at the
182  // polyPatch level
183 
184  const fvPatch& nonOverlapPatch = this->nonOverlapPatch();
185  const cyclicACMIFvPatch& nbrACMI = neighbPatch();
186  const fvPatch& nbrNonOverlapPatch = nbrACMI.nonOverlapPatch();
187 
188  resetPatchAreas(*this);
189  resetPatchAreas(nonOverlapPatch);
190  resetPatchAreas(nbrACMI);
191  resetPatchAreas(nbrNonOverlapPatch);
192 
193  // Scale the mesh flux
194 
195  const labelListList& newSrcAddr = AMI().srcAddress();
196  const labelListList& newTgtAddr = AMI().tgtAddress();
197 
198  const fvMesh& mesh = boundaryMesh().mesh();
199  surfaceScalarField& meshPhi = const_cast<fvMesh&>(mesh).setPhi();
200  surfaceScalarField::Boundary& meshPhiBf = meshPhi.boundaryFieldRef();
201 
202  // Note: phip and phiNonOverlapp will be different sizes if new faces
203  // have been added
204  scalarField& phip = meshPhiBf[cyclicACMIPolyPatch_.index()];
205  scalarField& phiNonOverlapp =
206  meshPhiBf[nonOverlapPatch.patch().index()];
207 
208  const auto& localFaces = cyclicACMIPolyPatch_.localFaces();
209  const auto& localPoints = cyclicACMIPolyPatch_.localPoints();
210 
211  forAll(phip, facei)
212  {
213  if (newSrcAddr[facei].empty())
214  {
215  // AMI patch with no connection to other coupled faces
216  phip[facei] = 0.0;
217  }
218  else
219  {
220  // Scale the mesh flux according to the area fraction
221  const face& fAMI = localFaces[facei];
222 
223  // Note: using raw point locations to calculate the geometric
224  // area - faces areas are currently scaled (decoupled from
225  // mesh points)
226  const scalar geomArea = fAMI.mag(localPoints);
227  phip[facei] *= magSf()[facei]/geomArea;
228  }
229  }
230 
231  forAll(phiNonOverlapp, facei)
232  {
233  const scalar w = 1.0 - cyclicACMIPolyPatch_.srcMask()[facei];
234  phiNonOverlapp[facei] *= w;
235  }
236 
237  scalarField& nbrPhip = meshPhiBf[nbrACMI.patch().index()];
238  scalarField& nbrPhiNonOverlapp =
239  meshPhiBf[nbrNonOverlapPatch.patch().index()];
240 
241  const auto& nbrLocalFaces = nbrACMI.patch().localFaces();
242  const auto& nbrLocalPoints = nbrACMI.patch().localPoints();
243 
244  forAll(nbrPhip, facei)
245  {
246  if (newTgtAddr[facei].empty())
247  {
248  nbrPhip[facei] = 0.0;
249  }
250  else
251  {
252  const face& fAMI = nbrLocalFaces[facei];
253 
254  // Note: using raw point locations to calculate the geometric
255  // area - faces areas are currently scaled (decoupled from
256  // mesh points)
257  const scalar geomArea = fAMI.mag(nbrLocalPoints);
258  nbrPhip[facei] *= nbrACMI.magSf()[facei]/geomArea;
259  }
260  }
261 
262  forAll(nbrPhiNonOverlapp, facei)
263  {
264  const scalar w = 1.0 - cyclicACMIPolyPatch_.tgtMask()[facei];
265  nbrPhiNonOverlapp[facei] *= w;
266  }
267 }
268 
269 
270 // ************************************************************************* //
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::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::face::mag
scalar mag(const UList< point > &p) const
Magnitude of face area.
Definition: faceI.H:113
Foam::coupledFvPatch::delta
virtual tmp< vectorField > delta() const =0
Return delta (P to N) vectors across coupled patch.
Definition: coupledFvPatch.C:46
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:415
Foam::cyclicACMIPolyPatch::tolerance
static scalar tolerance()
Overlap tolerance.
Definition: cyclicACMIPolyPatchI.H:65
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::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
Foam::cyclicACMIFvPatch::delta
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
Definition: cyclicACMIFvPatch.C:111
Foam::cyclicACMIFvPatch::coupled
virtual bool coupled() const
Definition: cyclicACMIFvPatch.C:105
Foam::cyclicACMIFvPatch::internalFieldTransfer
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
Definition: cyclicACMIFvPatch.C:165
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
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< vector >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::fvPatch::Cf
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:113
Foam::cyclicACMIFvPatch::movePoints
virtual void movePoints()
Correct patches after moving points.
Definition: cyclicACMIFvPatch.C:174
Foam::cyclicACMIFvPatch::nonOverlapPatch
virtual const fvPatch & nonOverlapPatch() const
Return non-overlapping fvPatch.
Definition: cyclicACMIFvPatch.H:131
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::PrimitivePatch::localFaces
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatch.C:297
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
DebugPout
#define DebugPout
Report an information message using Foam::Pout.
Definition: messageStream.H:365
Foam::PrimitivePatch::localPoints
const Field< point_type > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatch.C:339
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:66
Foam::polyPatch::faceCentres
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:313
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:46
Foam::List< labelList >
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::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:156
Foam::fvPatch::patch
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:157
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
transform.H
3D tensor transformation operations.
Foam::fvPatch::magSf
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:150
Foam::cyclicACMIFvPatch::makeWeights
void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: cyclicACMIFvPatch.C:57
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::patchIdentifier::name
const word & name() const
The patch name.
Definition: patchIdentifier.H:134
Foam::patchIdentifier::index
label index() const
The index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:158
Foam::fvPatch::makeWeights
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: fvPatch.C:164
Foam::cyclicACMIFvPatch
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI)
Definition: cyclicACMIFvPatch.H:53