filmFlux.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) 2019 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "filmFlux.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "surfaceInterpolate.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace functionObjects
39 {
40  defineTypeNameAndDebug(filmFlux, 0);
41  addToRunTimeSelectionTable(functionObject, filmFlux, dictionary);
42 }
43 }
44 
45 
47 Foam::functionObjects::filmFlux::filmModel()
48 {
49  return time_.objectRegistry::lookupObject<filmType>(filmName_);
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 (
57  const word& name,
58  const Time& runTime,
59  const dictionary& dict
60 )
61 :
63  filmName_("surfaceFilmProperties"),
64  resultName_(scopedName("filmFlux"))
65 {
66  read(dict);
67 }
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
73 {
75  {
76  dict.readIfPresent<word>("film", filmName_);
77  dict.readIfPresent<word>("result", resultName_);
78 
79  return true;
80  }
81 
82  return false;
83 }
84 
85 
87 {
88  const auto& model = filmModel();
89 
90  const fvMesh& filmMesh = model.regionMesh();
91 
92  auto* resultPtr = filmMesh.getObjectPtr<surfaceScalarField>(resultName_);
93 
94  if (!resultPtr)
95  {
96  resultPtr = new surfaceScalarField
97  (
98  IOobject
99  (
100  resultName_,
101  time_.timeName(),
102  filmMesh,
106  ),
107  filmMesh,
109  );
110 
111  filmMesh.objectRegistry::store(resultPtr);
112  }
113 
114  auto& result = *resultPtr;
115 
116  const surfaceScalarField& phi = model.phi();
117  const volScalarField& magSf = model.magSf();
118  const volScalarField::Internal& vol = filmMesh.V();
119 
120  volScalarField height
121  (
122  IOobject
123  (
124  "height",
125  time_.timeName(),
126  filmMesh,
130  ),
131  filmMesh,
134  );
135 
136  auto& heightc = height.ref();
137 
138  heightc = max(dimensionedScalar("eps", dimLength, ROOTVSMALL), vol/magSf());
139  height.correctBoundaryConditions();
141  result = phi/fvc::interpolate(height);
142 
143  return true;
144 }
145 
146 
148 {
149  const auto& filmMesh = filmModel().regionMesh();
150 
151  filmMesh.lookupObject<surfaceScalarField>(resultName_).write();
152 
153  return true;
154 }
155 
156 
157 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
Foam::surfaceFields.
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Kinematic form of single-cell layer surface film model.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
engineTime & runTime
Ignore writing from objectRegistry::writeObject()
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
virtual bool read(const dictionary &)
Read the field data.
Definition: filmFlux.C:65
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
virtual bool write()
Write the field.
Definition: filmFlux.C:140
filmFlux(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: filmFlux.C:49
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual bool execute()
Execute.
Definition: filmFlux.C:79
Base class for function objects, adding functionality to read/write state information (data required ...
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Request registration (bool: true)
Do not request registration (bool: false)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const Time & time_
Reference to the time database.