div.H
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) 2012-2016 OpenFOAM Foundation
9 Copyright (C) 2020-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::functionObjects::div
29
30Group
31 grpFieldFunctionObjects
32
33Description
34 Computes the divergence of an input field.
35
36 Operands:
37 \table
38 Operand | Type | Location
39 input | {surfaceScalar,volVector}Field | $FOAM_CASE/<time>/<inpField>
40 output file | - | -
41 output field | volScalarField | $FOAM_CASE/<time>/<outField>
42 \endtable
43
44Usage
45 Minimal example by using \c system/controlDict.functions:
46 \verbatim
47 div1
48 {
49 // Mandatory entries (unmodifiable)
50 type div;
51 libs (fieldFunctionObjects);
52
53 // Mandatory (inherited) entry (runtime modifiable)
54 field <field>;
55
56 // Inherited entries
57 ...
58 }
59 \endverbatim
60
61 where the entries mean:
62 \table
63 Property | Description | Type | Reqd | Deflt
64 type | Type name: div | word | yes | -
65 libs | Library name: fieldFunctionObjects | word | yes | -
66 field | Name of the operand field | word | yes | -
67 \endtable
68
69 The inherited entries are elaborated in:
70 - \link functionObject.H \endlink
71 - \link fieldExpression.H \endlink
72 - \link zoneSubSet.H \endlink
73
74 Minimal example by using the \c postProcess utility:
75 \verbatim
76 postProcess -func "div(<field>)"
77 \endverbatim
78
79Note
80 - To execute \c div function object on an input <field>, a numerical scheme
81 should be defined for \c div(<field>) in \c system/fvSchemes.divSchemes.
82 - Optionally the user can specify \c cellZones to create a sub-mesh for the
83 \c div calculation.
84
85See also
86 - Foam::functionObject
87 - Foam::functionObjects::fvMeshFunctionObject
88 - Foam::functionObjects::fieldExpression
89 - ExtendedCodeGuide::functionObjects::field::div
90
91SourceFiles
92 div.C
93 divTemplates.C
94
95\*---------------------------------------------------------------------------*/
96
97#ifndef functionObjects_div_H
98#define functionObjects_div_H
99
100#include "fieldExpression.H"
101
102// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
103
104namespace Foam
105{
106namespace functionObjects
107{
108
109/*---------------------------------------------------------------------------*\
110 Class div Declaration
111\*---------------------------------------------------------------------------*/
112
113class div
114:
115 public fieldExpression
116{
117 // Private Member Functions
118
119 //- Calculate the divergence of either a
120 //- volScalarField or a surfaceScalarField and register the result
121 template<class FieldType>
122 bool calcDiv();
123
124 //- Calculate the divergence field and return true if successful
125 virtual bool calc();
126
127 //- Helper function for writing submesh fields
128 template<class Type>
129 bool writeField();
130
131
132public:
133
134 //- Runtime type information
135 TypeName("div");
136
137
138 // Constructors
139
140 //- Construct from Time and dictionary
141 div
142 (
143 const word& name,
144 const Time& runTime,
145 const dictionary& dict
146 );
147
148 //- No copy construct
149 div(const div&) = delete;
150
151 //- No copy assignment
152 void operator=(const div&) = delete;
153
154
155 //- Destructor
156 virtual ~div() = default;
157
158
159 // Member Functions
160
161 //- Write the result field
162 virtual bool write();
163};
164
165
166// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167
168} // End namespace functionObjects
169} // End namespace Foam
170
171// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172
173#ifdef NoRepository
174 #include "divTemplates.C"
175#endif
176
177// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178
179#endif
180
181// ************************************************************************* //
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const word & name() const noexcept
Return the name of this functionObject.
Computes the divergence of an input field.
Definition: div.H:155
void operator=(const div &)=delete
No copy assignment.
TypeName("div")
Runtime type information.
virtual ~div()=default
Destructor.
div(const div &)=delete
No copy construct.
div(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: div.C:76
virtual bool write()
Write the result field.
Definition: div.C:58
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
A class for handling words, derived from Foam::string.
Definition: word.H:68
engineTime & runTime
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
Namespace for OpenFOAM.
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73