Dual.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-2016 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 "Dual.H"
27 #include "coupledPointPatchField.H"
28 
29 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const fvMesh& mesh
35 )
36 {
37  autoPtr<labelList> s(new labelList(2));
38  s()[0] = mesh.nCells();
39  s()[1] = mesh.nPoints();
40  return s;
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class Type>
48 (
49  const IOobject& io,
50  const dictionary& dict,
51  const fvMesh& mesh
52 )
53 :
54  AveragingMethod<Type>(io, dict, mesh, size(mesh)),
55  volumeCell_(mesh.V()),
56  volumeDual_(mesh.nPoints(), 0.0),
57  dataCell_(FieldField<Field, Type>::operator[](0)),
58  dataDual_(FieldField<Field, Type>::operator[](1)),
59  tetVertices_(3),
60  tetCoordinates_(4)
61 {
62  forAll(this->mesh_.C(), celli)
63  {
64  List<tetIndices> cellTets =
65  polyMeshTetDecomposition::cellTetIndices(this->mesh_, celli);
66  forAll(cellTets, tetI)
67  {
68  const tetIndices& tetIs = cellTets[tetI];
69  const face& f = this->mesh_.faces()[tetIs.face()];
70  const scalar v = tetIs.tet(this->mesh_).mag();
71  volumeDual_[f[tetIs.faceBasePt()]] += v;
72  volumeDual_[f[tetIs.facePtA()]] += v;
73  volumeDual_[f[tetIs.facePtB()]] += v;
74  }
75  }
76 
77  mesh.globalData().syncPointData
78  (
79  volumeDual_,
82  );
83 }
84 
85 
86 template<class Type>
88 (
89  const Dual<Type>& am
90 )
91 :
93  volumeCell_(am.volumeCell_),
94  volumeDual_(am.volumeDual_),
95  dataCell_(FieldField<Field, Type>::operator[](0)),
96  dataDual_(FieldField<Field, Type>::operator[](1)),
97  tetVertices_(am.tetVertices_),
98  tetCoordinates_(am.tetCoordinates_)
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
103 
104 template<class Type>
106 {}
107 
108 
109 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
110 
111 template<class Type>
113 (
114  const point position,
115  const tetIndices& tetIs
116 ) const
117 {
118  const face& f = this->mesh_.faces()[tetIs.face()];
119 
120  tetVertices_[0] = f[tetIs.faceBasePt()];
121  tetVertices_[1] = f[tetIs.facePtA()];
122  tetVertices_[2] = f[tetIs.facePtB()];
123 
124  tetIs.tet(this->mesh_).barycentric(position, tetCoordinates_);
125 
126  tetCoordinates_ = max(tetCoordinates_, scalar(0));
127 }
128 
129 
130 template<class Type>
132 {
133  this->mesh_.globalData().syncPointData
134  (
135  dataDual_,
136  plusEqOp<Type>(),
138  );
139 }
140 
141 
142 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
143 
144 template<class Type>
146 (
147  const point position,
148  const tetIndices& tetIs,
149  const Type& value
150 )
151 {
152  tetGeometry(position, tetIs);
153 
154  dataCell_[tetIs.cell()] +=
155  tetCoordinates_[0]*value
156  / (0.25*volumeCell_[tetIs.cell()]);
157 
158  for(label i = 0; i < 3; i ++)
159  {
160  dataDual_[tetVertices_[i]] +=
161  tetCoordinates_[i+1]*value
162  / (0.25*volumeDual_[tetVertices_[i]]);
163  }
164 }
165 
166 
167 template<class Type>
169 (
170  const point position,
171  const tetIndices& tetIs
172 ) const
173 {
174  tetGeometry(position, tetIs);
175 
176  return
177  tetCoordinates_[0]*dataCell_[tetIs.cell()]
178  + tetCoordinates_[1]*dataDual_[tetVertices_[0]]
179  + tetCoordinates_[2]*dataDual_[tetVertices_[1]]
180  + tetCoordinates_[3]*dataDual_[tetVertices_[2]];
181 }
182 
183 
184 template<class Type>
187 (
188  const point position,
189  const tetIndices& tetIs
190 ) const
191 {
192  tetGeometry(position, tetIs);
193 
194  const label celli(tetIs.cell());
195 
196  const tensor T
197  (
198  inv
199  (
200  tensor
201  (
202  this->mesh_.points()[tetVertices_[0]] - this->mesh_.C()[celli],
203  this->mesh_.points()[tetVertices_[1]] - this->mesh_.C()[celli],
204  this->mesh_.points()[tetVertices_[2]] - this->mesh_.C()[celli]
205  )
206  )
207  );
208 
209  const vector t( - T.T().x() - T.T().y() - T.T().z());
210 
211  const TypeGrad S
212  (
213  dataDual_[tetVertices_[0]],
214  dataDual_[tetVertices_[1]],
215  dataDual_[tetVertices_[2]]
216  );
217 
218  const Type s(dataCell_[celli]);
219 
220  return (T & S) + (t*s);
221 }
222 
223 
224 template<class Type>
226 {
227  syncDualData();
228 
230 }
231 
232 
233 template<class Type>
235 (
236  const AveragingMethod<scalar>& weight
237 )
238 {
239  syncDualData();
240 
242 }
243 
244 
245 template<class Type>
248 {
249  return tmp<Field<Type>>(dataCell_);
250 }
251 
252 
253 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
label cell() const
Return the cell.
Definition: tetIndicesI.H:30
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
label facePtA() const
Return face point A.
Definition: tetIndicesI.H:48
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
TypeGrad interpolateGrad(const point position, const tetIndices &tetIs) const
Interpolate gradient.
Definition: Dual.C:187
outerProduct< vector, Type >::type TypeGrad
Protected typedefs.
Dual(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Constructors.
Definition: Dual.C:48
Generic field type.
Definition: FieldField.H:51
label faceBasePt() const
Return the face base point.
Definition: tetIndicesI.H:42
void add(const point position, const tetIndices &tetIs, const Type &value)
Member Functions.
Definition: Dual.C:146
virtual ~Dual()
Destructor.
Definition: Dual.C:105
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dynamicFvMesh & mesh
label facePtB() const
Return face point B.
Definition: tetIndicesI.H:54
Default transformation behaviour.
Type interpolate(const point position, const tetIndices &tetIs) const
Interpolate.
Definition: Dual.C:169
Base class for lagrangian averaging methods.
Definition: MPPICParcel.H:55
List< label > labelList
A List of labels.
Definition: labelList.H:56
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
label face() const
Return the face.
Definition: tetIndicesI.H:36
Dual-mesh lagrangian averaging procedure.
Definition: Dual.H:64
scalar mag() const
Return volume.
Definition: tetrahedronI.H:171
labelList f(nPoints)
const volScalarField & T
AveragingMethod< Type >::TypeGrad TypeGrad
Public typedefs.
Definition: Dual.H:73
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
scalar barycentric(const point &pt, List< scalar > &bary) const
Calculate the barycentric coordinates of the given.
Definition: tetrahedronI.H:317
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
void average()
Calculate the average.
Definition: Dual.C:225
tmp< Field< Type > > primitiveField() const
Return an internal field of the average.
Definition: Dual.C:247
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
A class for managing temporary objects.
Definition: PtrList.H:54
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet from the.
Definition: tetIndicesI.H:66
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29