Loading...
Searching...
No Matches
Basic.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) 2013-2017 OpenFOAM Foundation
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\*---------------------------------------------------------------------------*/
28#include "Basic.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class Type>
35(
36 const IOobject& io,
37 const dictionary& dict,
38 const fvMesh& mesh
39)
40:
41 AveragingMethod<Type>(io, dict, mesh, labelList(1, mesh.nCells())),
42 data_(FieldField<Field, Type>::operator[](0)),
43 dataGrad_(mesh.nCells())
44{}
45
46
47template<class Type>
49(
50 const Basic<Type>& am
51)
52:
53 AveragingMethod<Type>(am),
54 data_(FieldField<Field, Type>::operator[](0)),
55 dataGrad_(am.dataGrad_)
56{}
57
58
59// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
60
61template<class Type>
63{}
64
65
66// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
67
68template<class Type>
70{
72 (
74 (
75 "BasicAverage::Data",
76 this->mesh_,
80 ),
81 this->mesh_,
84 );
85 tempData.primitiveFieldRef() = data_;
86 tempData.correctBoundaryConditions();
87 dataGrad_ = fvc::grad(tempData)->primitiveField();
88}
89
90
91// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92
93template<class Type>
95(
97 const tetIndices& tetIs,
98 const Type& value
100{
101 data_[tetIs.cell()] += value/this->mesh_.V()[tetIs.cell()];
102}
103
104
105template<class Type>
107(
109 const tetIndices& tetIs
110) const
112 return data_[tetIs.cell()];
113}
114
115
116template<class Type>
119(
121 const tetIndices& tetIs
122) const
124 return dataGrad_[tetIs.cell()];
125}
126
127
128template<class Type>
131{
132 return tmp<Field<Type>>(data_);
133}
134
135
136// ************************************************************************* //
Base class for lagrangian averaging methods.
const fvMesh & mesh_
The mesh on which the averaging is to be done.
AveragingMethod(const IOobject &io, const dictionary &dict, const fvMesh &mesh, const labelList &size)
Constructors.
virtual void updateGrad()
Protected member functions.
Basic lagrangian averaging procedure.
Definition Basic.H:64
AveragingMethod< Type >::TypeGrad TypeGrad
Public typedefs.
Definition Basic.H:74
Basic(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Constructors.
Definition Basic.C:28
tmp< Field< Type > > primitiveField() const
Return an internal field of the average.
Definition Basic.C:123
virtual ~Basic()
Destructor.
Definition Basic.C:55
Type interpolate(const barycentric &coordinates, const tetIndices &tetIs) const
Interpolate.
Definition Basic.C:100
void add(const barycentric &coordinates, const tetIndices &tetIs, const Type &value)
Member Functions.
Definition Basic.C:88
TypeGrad interpolateGrad(const barycentric &coordinates, const tetIndices &tetIs) const
Interpolate gradient.
Definition Basic.C:112
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
constexpr FieldField() noexcept
Construct null.
Definition FieldField.C:101
friend Ostream & operator(Ostream &, const FieldField< Field, Type > &)
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition tetIndices.H:79
label cell() const noexcept
Return the cell index.
Definition tetIndices.H:146
A class for managing temporary objects.
Definition tmp.H:75
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
const auto & io
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Definition List.H:62
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
dictionary dict