fvcSmooth.H
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 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10  Copyright (C) 2020 Henning Scheufler
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 InNamespace
29  Foam::fvc
30 
31 Description
32  Provides functions smooth spread and sweep which use the FaceCellWave
33  algorithm to smooth and redistribute the first field argument.
34 
35  smooth:
36  smooths the field by ensuring the values in neighbouring cells are
37  at least coeff* the cell value.
38 
39  spread:
40  redistributes the field by spreading the maximum value within the
41  region defined by the value (being between alphaMax and alphaMin)
42  and gradient of alpha (where the difference between the values in
43  neighbouring cells is larger than alphaDiff).
44 
45  sweep:
46  redistributes the field by sweeping the maximum value where the
47  gradient of alpha is large (where the difference between the values
48  in neighbouring cells is larger than alphaDiff) away from that
49  starting point of the sweep.
50 
51  spreadSource:
52  spread a source field (mDotIn) for two phase multiphase using
53  a laplacian operator and diffusivity D.
54  The spread source (mDotOut) is distributed from alpha1 < cutoff
55  to alpha1 > 1 - cutoff, and it is zero across the interface
56 
57 SourceFiles
58  fvcSmooth.C
59 
60 \*---------------------------------------------------------------------------*/
61 
62 #ifndef fvcSmooth_H
63 #define fvcSmooth_H
64 
65 #include "volFieldsFwd.H"
66 #include "dimensionedScalar.H"
67 
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70 namespace Foam
71 {
72 namespace fvc
73 {
74  void smooth
75  (
77  const scalar coeff
78  );
79 
80  void spread
81  (
83  const volScalarField& alpha,
84  const label nLayers,
85  const scalar alphaDiff = 0.2,
86  const scalar alphaMax = 0.99,
87  const scalar alphaMin = 0.01
88  );
89 
90  void sweep
91  (
93  const volScalarField& alpha,
94  const label nLayers,
95  const scalar alphaDiff = 0.2
96  );
97 
98  void spreadSource
99  (
100  volScalarField& mDotOut,
101  const volScalarField& mDotIn,
102  const volScalarField& alpha1,
103  const volScalarField& alpha2,
104  const dimensionedScalar& D,
105  const scalar cutoff
106  );
107 }
108 }
109 
110 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
111 
112 #endif
113 
114 // ************************************************************************* //
rDeltaTY field()
void spread(volScalarField &field, const volScalarField &alpha, const label nLayers, const scalar alphaDiff=0.2, const scalar alphaMax=0.99, const scalar alphaMin=0.01)
Definition: fvcSmooth.C:126
void spreadSource(volScalarField &mDotOut, const volScalarField &mDotIn, const volScalarField &alpha1, const volScalarField &alpha2, const dimensionedScalar &D, const scalar cutoff)
Definition: fvcSmooth.C:318
const volScalarField & alpha2
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:37
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
void sweep(volScalarField &field, const volScalarField &alpha, const label nLayers, const scalar alphaDiff=0.2)
Definition: fvcSmooth.C:224
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionedScalar & D
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.
const volScalarField & alpha1