ddt2Templates.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) 2016 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
28#include "Time.H"
29#include "dimensionedType.H"
30#include "fvcDdt.H"
31
32// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33
34template<class FieldType>
35int Foam::functionObjects::ddt2::apply(const word& inputName, int& state)
36{
37 // State: return 0 (not-processed), -1 (skip), +1 ok
38
39 // Already done, or not available
40 if (state || !foundObject<FieldType>(inputName))
41 {
42 return state;
43 }
44
45 const FieldType& input = lookupObject<FieldType>(inputName);
46
47 word outputName(resultName_);
48 outputName.replace("@@", inputName);
49
50 results_.set(outputName);
51
52 if (!foundObject<volScalarField>(outputName))
53 {
54 const dimensionSet dims =
55 (
56 mag_
57 ? mag(input.dimensions()/dimTime)
58 : magSqr(input.dimensions()/dimTime)
59 );
60
61 auto tddt2 = tmp<volScalarField>::New
62 (
63 IOobject
64 (
67 mesh_,
70 ),
71 mesh_,
73 );
74
75 store(outputName, tddt2);
76 }
77
78 volScalarField& output = lookupObjectRef<volScalarField>(outputName);
79
80 if (mag_)
81 {
83 }
84 else
85 {
87 }
88
89 // Could add additional statistics here
90 Log << type() << ' ' << name()
91 << " field " << outputName
92 << " average: " << gAverage(output) << endl;
93
94 state = +1;
95 return state;
96}
97
98
99// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
bool set(const Key &key)
Same as insert (no value to overwrite)
Definition: HashSet.H:197
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
virtual void apply()=0
Apply bins.
const word & name() const noexcept
Return the name of this functionObject.
virtual const word & type() const =0
Runtime type information.
const fvMesh & mesh_
Reference to the fvMesh.
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
const Time & time_
Reference to the time database.
string & replace(const std::string &s1, const std::string &s2, size_type pos=0)
Definition: string.C:108
word outputName("finiteArea-edges.obj")
Calculate the first temporal derivative.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static Istream & input(Istream &is, IntRange< T > &range)
Definition: IntRanges.C:55
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:66
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)