Loading...
Searching...
No Matches
fvcReconstructMag.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-2016 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 "fvcReconstruct.H"
29#include "fvMesh.H"
30#include "volFields.H"
31#include "surfaceFields.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40
41namespace fvc
42{
43
44// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45
47{
48 const fvMesh& mesh = ssf.mesh();
49
50 const labelUList& owner = mesh.owner();
51 const labelUList& neighbour = mesh.neighbour();
52
53 const volVectorField& C = mesh.C();
54 const surfaceVectorField& Cf = mesh.Cf();
55 const surfaceVectorField& Sf = mesh.Sf();
56 const surfaceScalarField& magSf = mesh.magSf();
57
58 tmp<volScalarField> treconField
59 (
61 (
63 (
64 "reconstruct("+ssf.name()+')',
65 ssf.instance(),
66 mesh,
69 ),
70 mesh,
73 )
74 );
75 scalarField& rf = treconField.ref();
76
77 forAll(owner, facei)
78 {
79 label own = owner[facei];
80 label nei = neighbour[facei];
81
82 rf[own] += (Sf[facei] & (Cf[facei] - C[own]))*ssf[facei]/magSf[facei];
83 rf[nei] -= (Sf[facei] & (Cf[facei] - C[nei]))*ssf[facei]/magSf[facei];
84 }
85
87
88 forAll(bsf, patchi)
89 {
90 const fvsPatchScalarField& psf = bsf[patchi];
91
92 const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
93 const vectorField& pCf = Cf.boundaryField()[patchi];
94 const vectorField& pSf = Sf.boundaryField()[patchi];
95 const scalarField& pMagSf = magSf.boundaryField()[patchi];
96
97 forAll(pOwner, pFacei)
98 {
99 label own = pOwner[pFacei];
100 rf[own] +=
101 (pSf[pFacei] & (pCf[pFacei] - C[own]))
102 *psf[pFacei]/pMagSf[pFacei];
103 }
104 }
105
106 rf /= mesh.V();
108 treconField.ref().correctBoundaryConditions();
109
110 return treconField;
111}
112
113
115{
117 (
118 fvc::reconstructMag(tssf())
119 );
120 tssf.clear();
121 return tvf;
122}
123
124
125// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
126
127} // End namespace fvc
128
129// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130
131} // End namespace Foam
132
133// ************************************************************************* //
Graphite solid properties.
Definition C.H:49
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
dynamicFvMesh & mesh
Reconstruct volField from a face flux field.
Namespace of functions to calculate explicit derivatives.
tmp< volScalarField > reconstructMag(const surfaceScalarField &)
Namespace for OpenFOAM.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvsPatchField< scalar > fvsPatchScalarField
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.