MPLIC.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2020-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "MPLIC.H"
27 #include "MPLICcell.H"
28 #include "volPointInterpolation.H"
29 #include "syncTools.H"
30 #include "slicedSurfaceFields.H"
31 #include "upwind.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
38 
41 }
42 
43 
44 // * * * * * * * * * * * * * * * Private Functions * * * * * * * * * * * * * //
45 
47 (
48  const label celli,
49  const scalarField& phi,
50  scalarField& alphaf,
51  boolList& correctedFaces,
52  const DynamicList<scalar>& cellAlphaf,
53  const fvMesh& mesh
54 ) const
55 {
56  // Face owners reference
57  const labelList& own = mesh.faceOwner();
58 
59  // The cell face labels
60  const labelList& cFaces = mesh.cells()[celli];
61 
62  // Fill the alphaf with surface interpolation in direction of the flow
63  forAll(cFaces, i)
64  {
65  const label facei = cFaces[i];
66  const scalar phiSigni = sign(phi[facei]);
67 
68  if
69  (
70  (own[facei] == celli && phiSigni == 1)
71  || (own[facei] != celli && phiSigni == -1)
72  )
73  {
74  alphaf[facei] = cellAlphaf[i];
75  correctedFaces[facei] = true;
76  }
77  }
78 }
79 
80 
82 (
83  const volScalarField& alpha,
84  const surfaceScalarField& phi,
85  scalarField& initAlphaf,
86  const bool unweighted,
87  const scalar tol,
88  const bool isMPLIC
89 ) const
90 {
91  // Finite volume mesh reference
92  const fvMesh& mesh = alpha.mesh();
93 
94  // Reference to primitive mesh
95  const primitiveMesh& primMesh = mesh;
96 
97  // Velocity field reference
98  const volVectorField& U
99  (
100  mesh.lookupObject<const volVectorField>
101  (
102  IOobject::groupName("U", phi.group())
103  )
104  );
105 
106  // Interpolate alpha from volume to the points of the mesh
107  const scalarField alphap
108  (
110  );
111 
112  // Interpolate U from cell centres to the points of the mesh
113  vectorField Up;
114 
115  if (!unweighted)
116  {
118  }
119 
120  // Flatten down phi flux field
121  const scalarField splicedPhi
122  (
124  (
125  IOobject
126  (
127  "splicedPhi",
128  mesh.time().name(),
129  mesh
130  ),
131  phi,
132  false
133  ).splice()
134  );
135 
136  scalarField alphaf(mesh.nFaces(), 0);
137 
138  // Mark which faces are corrected by MPLIC
139  boolList correctedFaces(mesh.nFaces(), false);
140 
141  // Construct class for cell cut
142  MPLICcell cutCell(unweighted, isMPLIC);
143 
144  // Loop through all the cells
145  forAll(mesh.cells(), celli)
146  {
147  if (alpha[celli] < (1 - tol) && alpha[celli] > tol)
148  {
149  // Store cell information
151  (
152  primMesh,
153  alphap,
154  Up,
155  alpha[celli],
156  U[celli],
157  celli
158  );
159 
160  // Volume ratio matching algorithm
161  if (cutCell.matchAlpha(cellInfo))
162  {
163  // Fill cutCell.alphaf() with face values from this cell
164  setCellAlphaf
165  (
166  celli,
167  splicedPhi,
168  alphaf,
169  correctedFaces,
170  cutCell.alphaf(),
171  mesh
172  );
173  }
174  }
175  }
176 
177  // Synchronise across the processor and cyclic patches
178  syncTools::syncFaceList(mesh, alphaf, plusEqOp<scalar>());
179  syncTools::syncFaceList(mesh, correctedFaces, plusEqOp<bool>());
180 
181  // Correct selected faces
182  forAll(correctedFaces, facei)
183  {
184  if (correctedFaces[facei])
185  {
186  initAlphaf[facei] = alphaf[facei];
187  }
188  }
189 
190  if (!mesh.conformal())
191  {
193  << "The " << type() << " scheme is not compatible with "
194  << "non-conformal meshes" << exit(FatalError);
195  }
196 
197  // Convert the alphaPhi spliced field into a surfaceScalarField
198  tmp<surfaceScalarField> tsplicedAlpha
199  (
201  (
202  "alphaf",
204  (
205  IOobject
206  (
207  "alphaf",
208  mesh.time().name(),
209  mesh
210  ),
211  mesh,
212  dimless,
213  initAlphaf,
214  false
215  ),
217  )
218  );
219  surfaceScalarField& splicedAlpha = tsplicedAlpha.ref();
220 
221  forAll(mesh.boundary(), patchi)
222  {
223  if (alpha.boundaryField()[patchi].fixesValue())
224  {
225  splicedAlpha.boundaryFieldRef()[patchi] =
227  }
228  }
229 
230  return tsplicedAlpha;
231 }
232 
233 
234 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
235 
237 (
238  const VolField<scalar>& vf
239 ) const
240 {
241  tmp<surfaceScalarField> tvff(upwind<scalar>(mesh(), phi_).interpolate(vf));
242 
243  scalarField splicedTvff
244  (
246  (
247  IOobject
248  (
249  "splicedTvff",
250  mesh().time().name(),
251  mesh()
252  ),
253  tvff,
254  false
255  ).splice()
256  );
257 
258  return surfaceAlpha(vf, phi_, splicedTvff, true, 1e-6);
259 }
260 
261 
262 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static volPointInterpolation & New(const word &name, const fvMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
static word groupName(Name name, const word &group)
Multicut Piecewise-Linear Interface Calculation (MPLIC) corrected scheme is a surface interpolation s...
Definition: MPLIC.H:86
virtual tmp< surfaceScalarField > interpolate(const VolField< scalar > &vf) const
Return the face-interpolate of the given cell field.
Definition: MPLIC.C:237
tmp< surfaceScalarField > surfaceAlpha(const volScalarField &alpha, const surfaceScalarField &phi, scalarField &splicedTvff, const bool unweighted, const scalar tol, const bool isMPLIC=true) const
Return alpha interpolation.
Definition: MPLIC.C:82
void setCellAlphaf(const label celli, const scalarField &phi, scalarField &alphaf, boolList &correctedFaces, const DynamicList< scalar > &cellAlphaf, const fvMesh &mesh) const
Set alphaPhi for the faces of the given cell.
Definition: MPLIC.C:47
Provides local cell addressing for geometry and data for MPLIC class.
Class performs geometric matching of volume fraction and calculates surface interpolation of volume f...
Definition: MPLICcell.H:56
bool matchAlpha(const MPLICcellStorage &cellInfo)
Match cell volume ratios.
Definition: MPLICcell.C:1027
const DynamicList< scalar > & alphaf() const
Return face volume fraction.
Definition: MPLICcellI.H:269
Specialisation of GeometricField which holds slices of given complete fields in a form that they act ...
tmp< Field< Type > > splice() const
Splice the sliced field and return the complete field.
Holds information regarding type of cell. Used in inside/outside determination in cellClassification.
Definition: cellInfo.H:67
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
bool conformal() const
Return whether the fvMesh is conformal with the polyMesh.
Definition: fvMesh.C:874
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:83
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1339
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
label nFaces() const
const cellList & cells() const
Abstract base class for surface interpolation schemes.
const fvMesh & mesh() const
Return mesh reference.
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronise values on all mesh faces.
Definition: syncTools.H:387
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
tmp< PointField< Type > > interpolate(const VolField< Type > &) const
Interpolate volField using inverse distance weighting.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const doubleScalar e
Definition: doubleScalar.H:106
dimensionedScalar sign(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
surfaceInterpolationScheme< scalar >::addMeshFluxConstructorToTable< MPLIC > addMPLICScalarMeshFluxConstructorToTable_
Definition: MPLIC.C:40
defineTypeNameAndDebug(combustionModel, 0)
error FatalError
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488