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 -------------------------------------------------------------------------------
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 \*---------------------------------------------------------------------------*/
27 
28 #include "Basic.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<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 
47 template<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 
61 template<class Type>
63 {}
64 
65 
66 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
67 
68 template<class Type>
70 {
72  (
73  IOobject
74  (
75  "BasicAverage::Data",
76  this->mesh_,
77  IOobject::NO_READ,
78  IOobject::NO_WRITE,
79  IOobject::NO_REGISTER
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 
93 template<class Type>
95 (
96  const barycentric& coordinates,
97  const tetIndices& tetIs,
98  const Type& value
99 )
100 {
101  data_[tetIs.cell()] += value/this->mesh_.V()[tetIs.cell()];
102 }
103 
104 
105 template<class Type>
107 (
108  const barycentric& coordinates,
109  const tetIndices& tetIs
110 ) const
111 {
112  return data_[tetIs.cell()];
113 }
114 
115 
116 template<class Type>
119 (
120  const barycentric& coordinates,
121  const tetIndices& tetIs
122 ) const
123 {
124  return dataGrad_[tetIs.cell()];
125 }
126 
127 
128 template<class Type>
131 {
132  return tmp<Field<Type>>(data_);
133 }
134 
135 
136 // ************************************************************************* //
void add(const barycentric &coordinates, const tetIndices &tetIs, const Type &value)
Member Functions.
Definition: Basic.C:88
dictionary dict
AveragingMethod< Type >::TypeGrad TypeGrad
Public typedefs.
Definition: Basic.H:74
const word zeroGradientType
A zeroGradient patch field type.
Basic(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Constructors.
Definition: Basic.C:28
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
tmp< Field< Type > > primitiveField() const
Return an internal field of the average.
Definition: Basic.C:123
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Generic dimensioned Type class.
const dimensionSet dimless
Dimensionless.
Type interpolate(const barycentric &coordinates, const tetIndices &tetIs) const
Interpolate.
Definition: Basic.C:100
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
TypeGrad interpolateGrad(const barycentric &coordinates, const tetIndices &tetIs) const
Interpolate gradient.
Definition: Basic.C:112
Basic lagrangian averaging procedure.
Definition: Basic.H:61
dynamicFvMesh & mesh
Generic templated field type.
Definition: Field.H:62
label cell() const noexcept
Return the cell index.
Definition: tetIndices.H:146
Base class for lagrangian averaging methods.
virtual ~Basic()
Destructor.
Definition: Basic.C:55
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:78
PtrList< coordinateSystem > coordinates(solidRegions.size())
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:51
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127