solverInfo.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) 2015-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "solverInfo.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace functionObjects
37 {
38  defineTypeNameAndDebug(solverInfo, 0);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
45 
47 {
49  {
50  return;
51  }
52 
53  if (writtenHeader_)
54  {
55  writeBreak(file());
56  }
57  else
58  {
59  writeHeader(os, "Solver information");
60  }
61 
62  writeCommented(os, "Time");
63 
64  for (const word& fieldName : fieldSet_.selectionNames())
65  {
66  writeFileHeader<scalar>(os, fieldName);
67  writeFileHeader<vector>(os, fieldName);
68  writeFileHeader<sphericalTensor>(os, fieldName);
69  writeFileHeader<symmTensor>(os, fieldName);
70  writeFileHeader<tensor>(os, fieldName);
71  }
72 
73  os << endl;
74 
75  writtenHeader_ = true;
76 }
77 
78 
80 (
81  const word& fieldName
82 )
83 {
84  if (!writeResidualFields_)
85  {
86  return;
87  }
88 
89  const word residualName
90  (
91  IOobject::scopedName("initialResidual", fieldName)
92  );
93 
94  auto* fieldPtr = mesh_.getObjectPtr<IOField<scalar>>(residualName);
95 
96  if (!fieldPtr)
97  {
98  auto* fieldPtr =
99  new IOField<scalar>
100  (
101  IOobject
102  (
103  residualName,
104  mesh_.time().timeName(),
105  mesh_,
109  ),
110  Field<scalar>(mesh_.nCells(), Zero)
111  );
112 
113  fieldPtr->store();
114 
115  residualFieldNames_.insert(residualName);
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
123 (
124  const word& name,
125  const Time& runTime,
126  const dictionary& dict
127 )
128 :
129  fvMeshFunctionObject(name, runTime, dict),
130  writeFile(obr_, name, typeName, dict),
131  fieldSet_(mesh_),
132  residualFieldNames_(),
133  writeResidualFields_(false),
134  initialised_(false)
135 {
136  read(dict);
137 }
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 {
145  {
146  initialised_ = false;
147 
148  fieldSet_.read(dict);
149 
150  writeResidualFields_ = dict.getOrDefault("writeResidualFields", false);
151 
152  residualFieldNames_.clear();
153 
154  return true;
155  }
156 
157  return false;
158 }
159 
160 
162 {
163  // Note: delaying initialisation until after first iteration so that
164  // we can find wildcard fields
165  if (!initialised_)
166  {
167  writeFileHeader(file());
168 
169  if (writeResidualFields_)
170  {
171  for (const word& fieldName : fieldSet_.selectionNames())
172  {
173  initialiseResidualField<scalar>(fieldName);
174  initialiseResidualField<vector>(fieldName);
175  initialiseResidualField<sphericalTensor>(fieldName);
176  initialiseResidualField<symmTensor>(fieldName);
177  initialiseResidualField<tensor>(fieldName);
178  }
179  }
180 
181  initialised_ = true;
182  }
183 
184  writeCurrentTime(file());
185 
186  for (const word& fieldName : fieldSet_.selectionNames())
187  {
188  updateSolverInfo<scalar>(fieldName);
189  updateSolverInfo<vector>(fieldName);
190  updateSolverInfo<sphericalTensor>(fieldName);
191  updateSolverInfo<symmTensor>(fieldName);
192  updateSolverInfo<tensor>(fieldName);
193  }
195  file() << endl;
196 
197  return true;
198 }
199 
200 
202 {
203  for (const word& residualName : residualFieldNames_)
204  {
205  const auto* residualPtr =
206  mesh_.findObject<IOField<scalar>>(residualName);
207 
208  if (residualPtr)
209  {
210  volScalarField residual
211  (
212  IOobject
213  (
214  residualName,
215  mesh_.time().timeName(),
216  mesh_,
220  ),
221  mesh_,
224  );
225 
226  residual.primitiveFieldRef() = *residualPtr;
227  residual.correctBoundaryConditions();
228 
229  residual.write();
230  }
231  }
232 
233  return true;
234 }
235 
236 
237 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:264
dictionary dict
Writes solver information for a list of user-specified fields.
Definition: solverInfo.H:164
defineTypeNameAndDebug(ObukhovLength, 0)
static bool initialised_(false)
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:339
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void writeBreak(Ostream &os) const
Write a break marker to the stream.
Definition: writeFile.C:362
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
wordList selectionNames() const
Return the current field selection, in sorted order.
bool writtenHeader_
Flag to identify whether the header has been written.
Definition: writeFile.H:157
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
Abstract base-class for Time/database function objects.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
virtual bool execute()
Execute solverInfo.
Definition: solverInfo.C:154
Macros for easy insertion into run-time selection tables.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
virtual bool updateSelection()
Update the selection using current contents of obr_.
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:313
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
void writeFileHeader(Ostream &os)
Output file header information.
Definition: solverInfo.C:39
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
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
solverFieldSelection fieldSet_
Names of operand fields.
Definition: solverInfo.H:176
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
virtual bool read(const dictionary &dict)
Read optional controls.
void createResidualField(const word &fieldName)
Create and store a residual field on the mesh database.
Definition: solverInfo.C:73
solverInfo(const solverInfo &)=delete
No copy construct.
virtual bool write()
Write solverInfo results.
Definition: solverInfo.C:194
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Request registration (bool: true)
A primitive field of type <T> with automated input and output.
Do not request registration (bool: false)
Namespace for OpenFOAM.
virtual bool read(const dictionary &)
Read solverInfo settings.
Definition: solverInfo.C:135
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127