Loading...
Searching...
No Matches
surfaceCourantNumber.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) 2024-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
29#include "faMesh.H"
30#include "fvMesh.H"
31#include "areaFields.H"
32#include "edgeFields.H"
33#include "facEdgeIntegrate.H"
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace functionObjects
41{
44 (
48 );
49}
50}
51
52
53// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54
55void Foam::functionObjects::surfaceCourantNumber::writeFileHeader(Ostream& os)
56{
57 writeHeader(os, "Surface Courant Number");
58
59 writeCommented(os, "Time");
60 writeTabbed(os, "min");
61 writeTabbed(os, "max");
62 writeTabbed(os, "mean");
63 os << endl;
65 writtenHeader_ = true;
66}
67
68
69// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70
72(
73 const word& name,
74 const Time& runTime,
75 const dictionary& dict
76)
77:
79 writeFile(mesh_, name, typeName, dict),
80 resultName_("surfaceCo"),
81 phisName_("phis"),
82 rhoName_("rho"),
83 faMeshPtr_(nullptr)
85 read(dict);
86}
87
88
89// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90
92{
94 {
95 return false;
96 }
97
98 dict.readIfPresent("result", resultName_);
99 dict.readIfPresent("phis", phisName_);
100 dict.readIfPresent("rho", rhoName_);
101
102 // Registry containing all finite-area meshes on the polyMesh
103 const auto* faRegistry = faMesh::registry(mesh_);
104
105 if (!faRegistry)
106 {
108 << "No finite-area object registry is available."
110 }
111
112 word areaName;
113 if (!dict.readIfPresent("area", areaName))
114 {
115 wordList available = faRegistry->sortedNames<faMesh>();
116 if (!available.empty())
117 {
118 areaName = available.front();
119 }
120 }
121
122 faMeshPtr_ = faRegistry->cfindObject<faMesh>(areaName);
123
124 if (!faMeshPtr_)
125 {
127 << "No finite-area mesh available."
129 }
130
131 return true;
132}
133
134
136{
137 const auto* phiPtr = faMeshPtr_->cfindObject<edgeScalarField>(phisName_);
138
139 if (!phiPtr)
140 {
142 << "No edge flux field is available. "
143 << "Name of provided edge flux field (phi): " << phisName_
144 << endl;
145
146 return false;
147 }
148
149 const auto& phi = *phiPtr;
150
151 tmp<areaScalarField::Internal> tCo =
152 (
153 (0.5*faMeshPtr_->time().deltaT())
154 * fac::edgeSum(Foam::mag(phi))().internalField()
155 / faMeshPtr_->S()
156 );
157
158 if (tCo().dimensions() == dimDensity)
159 {
160 tCo.ref() /= faMeshPtr_->lookupObject<areaScalarField>(rhoName_);
161 }
162
163
164 auto* resultPtr = faMeshPtr_->getObjectPtr<areaScalarField>(resultName_);
165
166 if (!resultPtr)
167 {
168 resultPtr = new areaScalarField
169 (
170 IOobject
171 (
172 resultName_,
173 faMeshPtr_->time().timeName(),
174 *faMeshPtr_,
175 IOobjectOption()
176 ),
177 *faMeshPtr_,
180 );
181 regIOobject::store(resultPtr);
182 }
183 auto& result = *resultPtr;
184
185 result.internalFieldRef() = tCo;
186 result.correctBoundaryConditions();
187
188
189 const scalarMinMax limits = gMinMax(result.primitiveField());
190 const scalar mean = gAverage(result.primitiveField());
191
192 Log << "Surface Courant number: "
193 << "mean: " << mean
194 << " max: " << limits.max()
195 << endl;
196
197 if (writeToFile())
198 {
199 if (!writtenHeader_) writeFileHeader(file());
200
201 writeCurrentTime(file());
202 file()
203 << token::TAB << limits.min()
204 << token::TAB << limits.max()
205 << token::TAB << mean
207 }
208
209 return true;
210}
211
212
214{
215 const auto* result = faMeshPtr_->cfindObject<areaScalarField>(resultName_);
216
217 if (!result)
218 {
219 return false;
220 }
221
222 Log << type() << " " << name() << " write: " << result->name() << endl;
223
224 result->write();
225
226 return true;
227}
228
229
230// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Internal & internalFieldRef(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
T & front()
Access first element of the list, position [0].
Definition UListI.H:239
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition faMesh.H:140
static const objectRegistry * registry(const polyMesh &pMesh)
Find the singleton parent registry (on the polyMesh) that contains all objects related to finite-area...
Definition faMesh.C:145
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Abstract base-class for Time/database function objects.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
Computes the surface Courant number field at finite-area face centres.
surfaceCourantNumber(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
virtual bool execute()
Calculate the Courant number field and return true if successful.
virtual bool write()
Write the result field.
virtual bool read(const dictionary &)
Read the surfaceCourantNumber data.
Base class for writing single files from the function objects.
Definition writeFile.H:113
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition writeFile.C:334
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition writeFile.C:344
virtual bool read(const dictionary &dict)
Read.
Definition writeFile.C:240
bool writtenHeader_
Flag to identify whether the header has been written.
Definition writeFile.H:157
virtual OFstream & file()
Return access to the file (if only 1).
Definition writeFile.C:270
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition writeFile.C:318
virtual void writeCurrentTime(Ostream &os) const
Write the current time to stream.
Definition writeFile.C:354
virtual bool writeToFile() const
Flag to allow writing to file.
Definition writeFile.C:286
bool store()
Register object with its registry and transfer ownership to the registry.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
@ TAB
Tab [isspace].
Definition token.H:142
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
auto limits
Definition setRDeltaT.H:186
engineTime & runTime
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
OBJstream os(runTime.globalPath()/outputName)
Edge integrate edgeField creating a areaField. Edge sum a edgeField creating a areaField.
auto & name
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< GeometricField< Type, faPatchField, areaMesh > > edgeSum(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
List< word > wordList
List of word.
Definition fileName.H:60
const dimensionSet dimless
Dimensionless.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
errorManip< error > abort(error &err)
Definition errorManip.H:139
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
dictionary dict