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