resolutionIndexModel.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) 2022-2024 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 "resolutionIndexModel.H"
29 #include "fvMesh.H"
30 #include "ListOps.H"
31 #include "turbulenceFields.H"
32 #include "turbulenceModel.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
41 }
42 
43 
44 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
45 
47 {
48  auto tV = volScalarField::New
49  (
50  "V",
52  mesh_,
53  dimVolume,
55  );
56 
57  tV.ref().primitiveFieldRef() = mesh_.V();
58  tV.ref().correctBoundaryConditions();
59 
60  return tV;
61 }
62 
63 
64 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 
67 (
68  const word& name,
69  const fvMesh& mesh,
70  const dictionary& dict
71 )
72 :
73  mesh_(mesh),
74  resultName_()
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
81 {
82  resultName_ = dict.getOrDefault("result", type());
83 
84  auto* indexPtr = mesh_.getObjectPtr<volScalarField>(resultName_);
85 
86  if (!indexPtr)
87  {
88  indexPtr = new volScalarField
89  (
90  IOobject
91  (
92  resultName_,
93  mesh_.time().timeName(),
94  mesh_.thisDb(),
98  ),
99  mesh_,
102  );
103 
104  regIOobject::store(indexPtr);
105  }
106 
107  return true;
108 }
109 
110 
111 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Various functions to operate on Lists.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
resolutionIndexModel(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
dynamicFvMesh & mesh
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
A base class for resolutionIndex models.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
tmp< volScalarField > V() const
Return cell volume field.
Reading is optional [identical to READ_IF_PRESENT].
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
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
virtual bool read(const dictionary &dict)
Read top-level dictionary.