Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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 "Dual.H"
30
31// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32
33template<class Type>
34Foam::labelList Foam::AveragingMethods::Dual<Type>::sizing(const fvMesh& mesh)
35{
36 labelList sizes(2);
37 sizes[0] = mesh.nCells();
38 sizes[1] = mesh.nPoints();
39 return sizes;
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45template<class Type>
47(
48 const IOobject& io,
49 const dictionary& dict,
50 const fvMesh& mesh
51)
52:
53 AveragingMethod<Type>(io, dict, mesh, sizing(mesh)),
54 volumeCell_(mesh.V()),
55 volumeDual_(mesh.nPoints(), Zero),
56 dataCell_(FieldField<Field, Type>::operator[](0)),
57 dataDual_(FieldField<Field, Type>::operator[](1))
58{
59 forAll(this->mesh_.C(), celli)
60 {
61 List<tetIndices> cellTets =
63 forAll(cellTets, tetI)
64 {
65 const tetIndices& tetIs = cellTets[tetI];
66 const triFace triIs = tetIs.faceTriIs(this->mesh_);
67 const scalar v = tetIs.tet(this->mesh_).mag();
68 volumeDual_[triIs[0]] += v;
69 volumeDual_[triIs[1]] += v;
70 volumeDual_[triIs[2]] += v;
71 }
72 }
73
74 mesh.globalData().syncPointData
75 (
76 volumeDual_,
79 );
80}
81
82
83template<class Type>
85(
86 const Dual<Type>& am
87)
88:
89 AveragingMethod<Type>(am),
90 volumeCell_(am.volumeCell_),
91 volumeDual_(am.volumeDual_),
92 dataCell_(FieldField<Field, Type>::operator[](0)),
93 dataDual_(FieldField<Field, Type>::operator[](1))
94{}
95
96
97// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
98
99template<class Type>
100void Foam::AveragingMethods::Dual<Type>::syncDualData()
101{
102 this->mesh_.globalData().syncPointData
103 (
104 dataDual_,
107 );
108}
109
110
111// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112
113template<class Type>
115(
117 const tetIndices& tetIs,
118 const Type& value
119)
120{
121 const triFace triIs(tetIs.faceTriIs(this->mesh_));
122
123 dataCell_[tetIs.cell()] +=
124 coordinates[0]*value
125 / (0.25*volumeCell_[tetIs.cell()]);
126
127 for (label i = 0; i < 3; ++i)
128 {
129 dataDual_[triIs[i]] +=
130 coordinates[i+1]*value
131 / (0.25*volumeDual_[triIs[i]]);
132 }
133}
134
135
136template<class Type>
138(
140 const tetIndices& tetIs
141) const
142{
143 const triFace triIs(tetIs.faceTriIs(this->mesh_));
144
145 return
146 coordinates[0]*dataCell_[tetIs.cell()]
147 + coordinates[1]*dataDual_[triIs[0]]
148 + coordinates[2]*dataDual_[triIs[1]]
149 + coordinates[3]*dataDual_[triIs[2]];
150}
151
152
153template<class Type>
156(
158 const tetIndices& tetIs
159) const
160{
161 const triFace triIs(tetIs.faceTriIs(this->mesh_));
162
163 const label celli(tetIs.cell());
164
165 const tensor T
166 (
167 inv
168 (
169 tensor
170 (
171 this->mesh_.points()[triIs[0]] - this->mesh_.C()[celli],
172 this->mesh_.points()[triIs[1]] - this->mesh_.C()[celli],
173 this->mesh_.points()[triIs[2]] - this->mesh_.C()[celli]
174 )
175 )
176 );
177
178 const vector t( - T.T().x() - T.T().y() - T.T().z());
179
180 const TypeGrad S
181 (
182 dataDual_[triIs[0]],
183 dataDual_[triIs[1]],
184 dataDual_[triIs[2]]
185 );
186
187 const Type s(dataCell_[celli]);
188
189 return (T & S) + (t*s);
190}
191
192
193template<class Type>
195{
196 syncDualData();
197
199}
200
201
202template<class Type>
204(
205 const AveragingMethod<scalar>& weight
206)
207{
208 syncDualData();
211}
212
213
214template<class Type>
217{
218 return tmp<Field<Type>>(dataCell_);
219}
220
221
222// ************************************************************************* //
Base class for lagrangian averaging methods.
const fvMesh & mesh_
The mesh on which the averaging is to be done.
AveragingMethod(const IOobject &io, const dictionary &dict, const fvMesh &mesh, const labelList &size)
Constructors.
virtual void average()
Calculate the average.
Dual-mesh lagrangian averaging procedure.
Definition Dual.H:68
AveragingMethod< Type >::TypeGrad TypeGrad
Public typedefs.
Definition Dual.H:78
tmp< Field< Type > > primitiveField() const
Return an internal field of the average.
Definition Dual.C:209
Dual(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Constructors.
Definition Dual.C:40
Type interpolate(const barycentric &coordinates, const tetIndices &tetIs) const
Interpolate.
Definition Dual.C:131
void add(const barycentric &coordinates, const tetIndices &tetIs, const Type &value)
Member Functions.
Definition Dual.C:108
TypeGrad interpolateGrad(const barycentric &coordinates, const tetIndices &tetIs) const
Interpolate gradient.
Definition Dual.C:149
void average()
Calculate the average.
Definition Dual.C:187
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
tmp< FieldField< Field, Type > > T() const
Return the field transpose (only defined for second rank tensors).
Definition FieldField.C:366
constexpr FieldField() noexcept
Construct null.
Definition FieldField.C:101
friend Ostream & operator(Ostream &, const FieldField< Field, Type > &)
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const volVectorField & C() const
Return cell centres as volVectorField.
Default transformation behaviour.
static List< tetIndices > cellTetIndices(const polyMesh &mesh, label cI)
Return the tet decomposition of the given cell, see.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition tetIndices.H:79
triFace faceTriIs(const polyMesh &mesh, const bool warn=true) const
The indices corresponding to the tri on the face for this tet. The normal of the tri points out of th...
Definition tetIndicesI.H:48
tetPointRef tet(const polyMesh &mesh) const
The tet geometry for this tet, where point0 is the cell centre.
label cell() const noexcept
Return the cell index.
Definition tetIndices.H:146
scalar mag() const
Return volume.
A class for managing temporary objects.
Definition tmp.H:75
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition triFace.H:68
const volScalarField & T
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
const auto & io
label nPoints
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
List< label > labelList
A List of labels.
Definition List.H:62
Tensor< scalar > tensor
Definition symmTensor.H:57
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
Vector< scalar > vector
Definition vector.H:57
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299