Basic.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 -------------------------------------------------------------------------------
10 License
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 Class
27  Foam::AveragingMethods::Basic
28 
29 Group
30  grpLagrangianIntermediateMPPICAveragingMethods
31 
32 Description
33  Basic lagrangian averaging procedure.
34 
35  This is a cell-volume based average. Point values are summed over the
36  computational cells and the result is divided by the cell volume.
37 
38  Interpolation is done assuming a constant value over each cells. Cell
39  gradients are calculated by the default fvc::grad scheme, and are also
40  assumed constant when interpolated.
41 
42 SourceFiles
43  Basic.C
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #ifndef Basic_H
48 #define Basic_H
49 
50 #include "AveragingMethod.H"
51 #include "pointMesh.H"
52 #include "tetIndices.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 namespace AveragingMethods
59 {
60 
61 /*---------------------------------------------------------------------------*\
62  Class Basic Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 template<class Type>
66 class Basic
67 :
68  public AveragingMethod<Type>
69 {
70 public:
71 
72  //- Public typedefs
73 
74  //- Gradient type
76 
77 
78 private:
79 
80  //- Private data
81 
82  //- Cell average field
83  Field<Type>& data_;
84 
85  //- Gradient field
86  Field<TypeGrad> dataGrad_;
87 
88 
89  //- Private member functions
90 
91  //- Re-calculate gradient
92  virtual void updateGrad();
93 
94 
95 public:
96 
97  //- Runtime type information
98  TypeName("basic");
99 
100 
101  //- Constructors
102 
103  //- Construct from components
104  Basic
105  (
106  const IOobject& io,
107  const dictionary& dict,
108  const fvMesh &mesh
109  );
110 
111  //- Construct a copy
112  Basic(const Basic<Type>& am);
113 
114  //- Construct and return a clone
115  virtual autoPtr<AveragingMethod<Type>> clone() const
116  {
118  (
119  new Basic<Type>(*this)
120  );
121  }
122 
123 
124  //- Destructor
125  virtual ~Basic();
126 
127 
128  //- Member Functions
129 
130  //- Add point value to interpolation
131  void add
132  (
133  const barycentric& coordinates,
134  const tetIndices& tetIs,
135  const Type& value
136  );
137 
138  //- Interpolate
139  Type interpolate
140  (
141  const barycentric& coordinates,
142  const tetIndices& tetIs
143  ) const;
144 
145  //- Interpolate gradient
147  (
148  const barycentric& coordinates,
149  const tetIndices& tetIs
150  ) const;
151 
152  //- Return an internal field of the average
154 
155  //- Return an internal field of the gradient
157 };
158 
159 
160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161 
162 } // End namespace AveragingMethods
163 } // End namespace Foam
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
167 #ifdef NoRepository
168  #include "Basic.C"
169 #endif
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 #endif
174 
175 // ************************************************************************* //
Foam::AveragingMethod
Base class for lagrangian averaging methods.
Definition: KinematicParcel.H:69
Foam::AveragingMethods::Basic::internalFieldGrad
tmp< Field< TypeGrad > > internalFieldGrad() const
Return an internal field of the gradient.
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::AveragingMethods::Basic::Basic
Basic(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Constructors.
Definition: Basic.C:35
AveragingMethod.H
Foam::AveragingMethods::Basic::add
void add(const barycentric &coordinates, const tetIndices &tetIs, const Type &value)
Member Functions.
Definition: Basic.C:95
Foam::Field
Generic templated field type.
Definition: Field.H:63
coordinates
PtrList< coordinateSystem > coordinates(solidRegions.size())
Foam::AveragingMethods::Basic::primitiveField
tmp< Field< Type > > primitiveField() const
Return an internal field of the average.
Definition: Basic.C:130
Foam::Barycentric< scalar >
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::AveragingMethods::Basic::interpolateGrad
TypeGrad interpolateGrad(const barycentric &coordinates, const tetIndices &tetIs) const
Interpolate gradient.
Definition: Basic.C:119
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::AveragingMethods::Basic::~Basic
virtual ~Basic()
Destructor.
Definition: Basic.C:62
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:83
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Basic.C
Foam::AveragingMethods::Basic::interpolate
Type interpolate(const barycentric &coordinates, const tetIndices &tetIs) const
Interpolate.
Definition: Basic.C:107
Foam::AveragingMethods::Basic
Basic lagrangian averaging procedure.
Definition: Basic.H:65
Foam::AveragingMethods::Basic::TypeName
TypeName("basic")
Runtime type information.
Foam::AveragingMethods::Basic::clone
virtual autoPtr< AveragingMethod< Type > > clone() const
Construct and return a clone.
Definition: Basic.H:114
Foam::AveragingMethods::Basic::TypeGrad
AveragingMethod< Type >::TypeGrad TypeGrad
Public typedefs.
Definition: Basic.H:74
tetIndices.H
pointMesh.H
Foam::AveragingMethod::TypeGrad
outerProduct< vector, Type >::type TypeGrad
Protected typedefs.
Definition: AveragingMethod.H:68