weightedFlux.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) 2019 Norbert Weber, HZDR
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 "weightedFlux.H"
29#include "demandDrivenData.H"
30
31// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
32
33template<class Type>
35{
38}
39
40
41// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
42
43template<class Type>
45{
46 clearOut();
47}
48
49
50// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51
52template<class Type>
54{
55 const fvMesh& mesh = this->mesh();
56
57 oDelta_ = new surfaceScalarField
58 (
60 (
61 "oDelta",
62 mesh.pointsInstance(),
63 mesh
64 ),
65 mesh,
66 dimLength
67 );
68 auto& oDelta = *oDelta_;
69
70 nDelta_ = new surfaceScalarField
71 (
73 (
74 "nDelta",
75 mesh.pointsInstance(),
76 mesh
77 ),
78 mesh,
79 dimLength
80 );
81 auto& nDelta = *nDelta_;
82
83 const labelUList& owner = mesh.owner();
84 const labelUList& neighbour = mesh.neighbour();
85
86 const surfaceVectorField n(mesh.Sf()/mesh.magSf());
87
88 const vectorField& C = mesh.cellCentres();
89 const vectorField& Cf = mesh.faceCentres();
90
91 // all distances are NORMAL to the face,
92 // as in the weighting factors in surfaceInterpolation.C
93 forAll(owner, facei)
94 {
95 oDelta[facei] =
96 mag(n[facei] & (C[owner[facei]] - Cf[facei]));
97 nDelta[facei] =
98 mag(n[facei] & (C[neighbour[facei]] - Cf[facei]));
99 }
100
101 const fvPatchList& patches = mesh.boundary();
102
103 forAll(patches, patchi)
104 {
105 const fvPatch& currPatch = mesh.boundary()[patchi];
106
107 // Patch normal vector
108 const vectorField nPatch(currPatch.Sf()/currPatch.magSf());
109
110 // Processor patch
111 if (currPatch.coupled())
112 {
113 const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
114 const vectorField& pCf = mesh.Cf().boundaryField()[patchi];
115
116 forAll(pOwner, facei)
117 {
118 const label own = pOwner[facei];
119
120 // All distances are NORMAL to the face
121 oDelta.boundaryFieldRef()[patchi][facei] =
122 mag(nPatch[facei] & (pCf[facei] - C[own]));
123 }
124
125 // Weight = delta_neighbour / delta in ORTHOGONAL direction,
126 nDelta.boundaryFieldRef()[patchi] =
127 currPatch.weights()*oDelta.boundaryFieldRef()[patchi]
128 /(scalar(1) - currPatch.weights());
129 }
130 else
131 {
132 const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
133 const vectorField& pCf = mesh.Cf().boundaryField()[patchi];
134
135 forAll(pOwner, facei)
136 {
137 const label own = pOwner[facei];
138
139 // All distances are NORMAL to the face!
140 oDelta.boundaryFieldRef()[patchi][facei] =
141 mag(nPatch[facei] & (pCf[facei] - C[own]));
142
143 nDelta.boundaryFieldRef()[patchi][facei] =
144 mag(nPatch[facei] & (pCf[facei] - C[own]));
145 }
146 }
147 }
148}
149
150
151template<class Type>
154(
156) const
157{
158 const fvMesh& mesh = vf.mesh();
161
163 (
165 (
166 "weightedFlux::interpolate(" + vf.name() + ')',
167 mesh.time().timeName(),
168 mesh
169 ),
170 mesh,
171 vf.dimensions()
172 );
173 auto& result = tresult.ref();
174
175 const labelUList& owner = mesh.owner();
176 const labelUList& neighbour = mesh.neighbour();
177
178 forAll(result, facei)
179 {
180 const scalar sigmaDeltaO = sigma_[owner[facei]]/oDelta[facei];
181 const scalar sigmaDeltaN = sigma_[neighbour[facei]]/nDelta[facei];
182
183 result[facei] =
184 (vf[owner[facei]]*sigmaDeltaO + vf[neighbour[facei]]*sigmaDeltaN)
185 /(sigmaDeltaO + sigmaDeltaN);
186 }
187
188 // Interpolate patches
189 auto& bfld = result.boundaryFieldRef();
190
191 forAll(bfld, patchi)
192 {
193 fvsPatchField<Type>& pfld = bfld[patchi];
194
195 // If not coupled - simply copy the boundary values of the field
196 if (!pfld.coupled())
197 {
198 pfld = vf.boundaryField()[patchi];
199 }
200 else
201 {
202 // e.g. processor patches have to calculated separately
203
204 const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
205 scalarField sigmaN
206 (
207 sigma_.boundaryField()[patchi].patchNeighbourField()
208 );
209
210 Field<Type> vfO(vf.boundaryField()[patchi].patchInternalField());
211 Field<Type> vfN(vf.boundaryField()[patchi].patchNeighbourField());
212
213 forAll(pOwner, facei)
214 {
215 const label own = pOwner[facei];
216
217 const scalar sigmaDeltaO =
218 sigma_[own]/oDelta.boundaryField()[patchi][facei];
219
220 const scalar sigmaDeltaN =
221 sigmaN[facei]/nDelta.boundaryField()[patchi][facei];
222
223 pfld[facei] =
224 (vfO[facei]*sigmaDeltaO + vfN[facei]*sigmaDeltaN)
225 /(sigmaDeltaO + sigmaDeltaN);
226 }
227 }
228 }
229
230 return tresult;
231}
232
233
234// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235
236namespace Foam
237{
239}
240
241
242// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
label n
Graphite solid properties.
Definition: C.H:53
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
virtual bool coupled() const
Return true if this patch field is coupled.
bool interpolate() const noexcept
Same as isPointData()
A class for managing temporary objects.
Definition: tmp.H:65
Weighted flux interpolation scheme class.
Definition: weightedFlux.H:89
virtual ~weightedFlux()
Destructor.
Definition: weightedFlux.C:44
void clearOut()
Clear all fields.
Definition: weightedFlux.C:34
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
Template functions to aid in the implementation of demand driven data.
Namespace for OpenFOAM.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
void deleteDemandDrivenData(DataPtr &dataPtr)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define makeSurfaceInterpolationScheme(SS)