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 | Copyright (C) 2013 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"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class Type>
32 (
33  const IOobject& io,
34  const dictionary& dict,
35  const fvMesh& mesh
36 )
37 :
38  AveragingMethod<Type>(io, dict, mesh, labelList(1, mesh.nCells())),
40  dataGrad_(mesh.nCells())
41 {}
42 
43 
44 template<class Type>
46 (
47  const Basic<Type>& am
48 )
49 :
52  dataGrad_(am.dataGrad_)
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
57 
58 template<class Type>
60 {}
61 
62 
63 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
64 
65 template<class Type>
67 {
69  (
70  IOobject
71  (
72  "BasicAverage::Data",
73  this->mesh_,
74  IOobject::NO_READ,
75  IOobject::NO_WRITE,
76  false
77  ),
78  this->mesh_,
81  );
82  tempData.internalField() = data_;
83  tempData.correctBoundaryConditions();
84  dataGrad_ = fvc::grad(tempData)->internalField();
85 }
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
90 template<class Type>
92 (
93  const point position,
94  const tetIndices& tetIs,
95  const Type& value
96 )
97 {
98  data_[tetIs.cell()] += value/this->mesh_.V()[tetIs.cell()];
99 }
100 
101 
102 template<class Type>
104 (
105  const point position,
106  const tetIndices& tetIs
107 ) const
108 {
109  return data_[tetIs.cell()];
110 }
111 
112 
113 template<class Type>
116 (
117  const point position,
118  const tetIndices& tetIs
119 ) const
120 {
121  return dataGrad_[tetIs.cell()];
122 }
123 
124 
125 template<class Type>
128 {
129  return tmp<Field<Type> >(data_);
130 }
131 
132 
133 // ************************************************************************* //
label cell() const
Return the cell.
Definition: tetIndicesI.H:30
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Generic field type.
Definition: FieldField.H:51
InternalField & internalField()
Return internal field.
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual ~Basic()
Destructor.
Definition: Basic.C:59
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
TypeGrad interpolateGrad(const point position, const tetIndices &tetIs) const
Interpolate gradient.
Definition: Basic.C:116
dictionary dict
Generic dimensioned Type class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Basic lagrangian averaging procedure.
Definition: Basic.H:60
Type interpolate(const point position, const tetIndices &tetIs) const
Interpolate.
Definition: Basic.C:104
Basic(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Constructors.
Definition: Basic.C:32
AveragingMethod< Type >::TypeGrad TypeGrad
Public typedefs.
Definition: Basic.H:69
Generic GeometricField class.
Traits class for primitives.
Definition: pTraits.H:50
void correctBoundaryConditions()
Correct boundary field.
List< label > labelList
A List of labels.
Definition: labelList.H:56
tmp< Field< Type > > internalField() const
Return an internal field of the average.
Definition: Basic.C:127
Base class for lagrangian averaging methods.
Definition: MPPICParcel.H:55
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
void add(const point position, const tetIndices &tetIs, const Type &value)
Member Functions.
Definition: Basic.C:92
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:73
A class for managing temporary objects.
Definition: PtrList.H:118