Basic.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "Basic.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const IOobject& io,
35  const dictionary& dict,
36  const fvMesh& mesh
37 )
38 :
39  AveragingMethod<Type>(io, dict, mesh, labelList(1, mesh.nCells())),
41  dataGrad_(mesh.nCells())
42 {}
43 
44 
45 template<class Type>
47 (
48  const Basic<Type>& am
49 )
50 :
53  dataGrad_(am.dataGrad_)
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58 
59 template<class Type>
61 {}
62 
63 
64 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
65 
66 template<class Type>
68 {
70  (
71  IOobject
72  (
73  "BasicAverage::Data",
74  this->mesh_,
75  IOobject::NO_READ,
76  IOobject::NO_WRITE,
77  false
78  ),
79  this->mesh_,
80  dimensioned<Type>("zero", dimless, Zero),
82  );
83  tempData.primitiveFieldRef() = data_;
84  tempData.correctBoundaryConditions();
85  dataGrad_ = fvc::grad(tempData)->primitiveField();
86 }
87 
88 
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 
91 template<class Type>
93 (
94  const barycentric& coordinates,
95  const tetIndices& tetIs,
96  const Type& value
97 )
98 {
99  data_[tetIs.cell()] += value/this->mesh_.V()[tetIs.cell()];
100 }
101 
102 
103 template<class Type>
105 (
106  const barycentric& coordinates,
107  const tetIndices& tetIs
108 ) const
109 {
110  return data_[tetIs.cell()];
111 }
112 
113 
114 template<class Type>
117 (
118  const barycentric& coordinates,
119  const tetIndices& tetIs
120 ) const
121 {
122  return dataGrad_[tetIs.cell()];
123 }
124 
125 
126 template<class Type>
129 {
130  return tmp<Field<Type>>(data_);
131 }
132 
133 
134 // ************************************************************************* //
void add(const barycentric &coordinates, const tetIndices &tetIs, const Type &value)
Member Functions.
Definition: Basic.C:93
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
AveragingMethod< Type >::TypeGrad TypeGrad
Gradient type.
Definition: Basic.H:69
Basic(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Constructors.
Definition: Basic.C:33
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
tmp< Field< Type > > primitiveField() const
Return an internal field of the average.
Definition: Basic.C:128
Generic GeometricField class.
const dimensionSet dimless
Generic dimensioned Type class.
Type interpolate(const barycentric &coordinates, const tetIndices &tetIs) const
Interpolate.
Definition: Basic.C:105
Generic field type.
Definition: FieldField.H:51
TypeGrad interpolateGrad(const barycentric &coordinates, const tetIndices &tetIs) const
Interpolate gradient.
Definition: Basic.C:117
Basic lagrangian averaging procedure.
Definition: Basic.H:60
dynamicFvMesh & mesh
label cell() const
Return the cell.
Definition: tetIndicesI.H:28
Base class for lagrangian averaging methods.
Definition: MPPICParcel.H:55
virtual ~Basic()
Destructor.
Definition: Basic.C:60
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void correctBoundaryConditions()
Correct boundary field.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92