Loading...
Searching...
No Matches
Dual.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) 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
26Class
27 Foam::AveragingMethods::Dual
28
29Group
30 grpLagrangianIntermediateMPPICAveragingMethods
31
32Description
33 Dual-mesh lagrangian averaging procedure.
34
35 Point values are summed using the tetrahedral decomposition of the
36 computational cells. Summation is done in the cells, and also in the
37 terahedrons surrounding each point. The latter forms a type of dual mesh.
38 The interpolation is weighted by proximity to the cell centre or point, as
39 calculated by the barycentric coordinate within the tethrahedron.
40
41 Values are interpolated linearly across the tethrahedron. Gradients are
42 calculated directly from the point values using a first order finite
43 element basis. The computed gradient is assumed constant over the
44 tethrahedron.
45
46SourceFiles
47 Dual.C
48
49\*---------------------------------------------------------------------------*/
50
51#ifndef Dual_H
52#define Dual_H
53
54#include "AveragingMethod.H"
55#include "pointMesh.H"
56#include "tetIndices.H"
57
58// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59
60namespace Foam
61{
62namespace AveragingMethods
63{
64
65/*---------------------------------------------------------------------------*\
66 Class Dual Declaration
67\*---------------------------------------------------------------------------*/
68
69template<class Type>
70class Dual
71:
72 public AveragingMethod<Type>
73{
74public:
75
76 //- Public typedefs
77
78 //- Gradient type
80
81
82private:
83
84 //- Private data
85
86 //- Volume of the cell-centered regions
87 const Field<scalar>& volumeCell_;
88
89 //- Volume of the point-centered regions
90 Field<scalar> volumeDual_;
91
92 //- Data on the cells
93 Field<Type>& dataCell_;
94
95 //- Data on the points
96 Field<Type>& dataDual_;
97
98
99 //- Private Member Functions
100
101 //- The sizing for the FieldField parts
102 static labelList sizing(const fvMesh& mesh);
103
104
105 //- Sync point data over processor boundaries
106 void syncDualData();
107
108
109public:
110
111 //- Runtime type information
112 TypeName("dual");
113
114
115 //- Constructors
116
117 //- Construct from components
118 Dual
119 (
120 const IOobject& io,
121 const dictionary& dict,
122 const fvMesh& mesh
123 );
124
125 //- Construct a copy
126 Dual(const Dual<Type>& am);
127
128 //- Construct and return a clone
130 {
132 (
133 new Dual<Type>(*this)
134 );
135 }
136
137
138 //- Destructor
139 virtual ~Dual() = default;
140
141
142 //- Member Functions
143
144 //- Add point value to interpolation
145 void add
146 (
148 const tetIndices& tetIs,
149 const Type& value
150 );
151
152 //- Interpolate
153 Type interpolate
156 const tetIndices& tetIs
157 ) const;
158
159 //- Interpolate gradient
161 (
163 const tetIndices& tetIs
164 ) const;
165
166 //- Calculate the average
167 void average();
168 void average(const AveragingMethod<scalar>& weight);
169
170 //- Return an internal field of the average
172
173 //- Return an internal field of the gradient
175};
176
177
178// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179
180} // End namespace AveragingMethods
181} // End namespace Foam
182
183// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184
185#ifdef NoRepository
186 #include "Dual.C"
187#endif
188
189// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190
191#endif
192
193// ************************************************************************* //
AveragingMethod(const IOobject &io, const dictionary &dict, const fvMesh &mesh, const labelList &size)
Constructors.
outerProduct< vector, Type >::type TypeGrad
Protected typedefs.
AveragingMethod< Type >::TypeGrad TypeGrad
Public typedefs.
Definition Dual.H:78
tmp< Field< TypeGrad > > internalFieldGrad() const
Return an internal field of the gradient.
virtual autoPtr< AveragingMethod< Type > > clone() const
Construct and return a clone.
Definition Dual.H:154
tmp< Field< Type > > primitiveField() const
Return an internal field of the average.
Definition Dual.C:209
Dual(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Constructors.
Definition Dual.C:40
Type interpolate(const barycentric &coordinates, const tetIndices &tetIs) const
Interpolate.
Definition Dual.C:131
void add(const barycentric &coordinates, const tetIndices &tetIs, const Type &value)
Member Functions.
Definition Dual.C:108
TypeGrad interpolateGrad(const barycentric &coordinates, const tetIndices &tetIs) const
Interpolate gradient.
Definition Dual.C:149
void average()
Calculate the average.
Definition Dual.C:187
virtual ~Dual()=default
Destructor.
TypeName("dual")
Runtime type information.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition tetIndices.H:79
A class for managing temporary objects.
Definition tmp.H:75
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
const auto & io
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68