snappyLayerDriverTemplates.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) 2011-2012 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 "snappyLayerDriver.H"
29#include "syncTools.H"
30
31// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32
33template<class Type>
34void Foam::snappyLayerDriver::averageNeighbours
35(
36 const polyMesh& mesh,
37 const bitSet& isMasterEdge,
38 const labelList& meshEdges,
39 const labelList& meshPoints,
40 const edgeList& edges,
41 const scalarField& invSumWeight,
42 const Field<Type>& data,
43 Field<Type>& average
44)
45{
46 const pointField& pts = mesh.points();
47
48 average = Type(Zero);
49
50 forAll(edges, edgeI)
51 {
52 if (isMasterEdge.test(meshEdges[edgeI]))
53 {
54 const edge& e = edges[edgeI];
55 //scalar eWeight = edgeWeights[edgeI];
56 //scalar eWeight = 1.0;
57 scalar eMag = max
58 (
59 VSMALL,
60 mag
61 (
62 pts[meshPoints[e[1]]]
63 - pts[meshPoints[e[0]]]
64 )
65 );
66 scalar eWeight = 1.0/eMag;
67
68 label v0 = e[0];
69 label v1 = e[1];
70
71 average[v0] += eWeight*data[v1];
72 average[v1] += eWeight*data[v0];
73 }
74 }
75
77 (
78 mesh,
79 meshPoints,
80 average,
81 plusEqOp<Type>(),
82 Type(Zero) // null value
83 );
84
85 average *= invSumWeight;
86}
87
88
89// ************************************************************************* //
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
dynamicFvMesh & mesh
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333