fvcReconstructMag.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fvcReconstruct.H"
27 #include "fvMesh.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace fvc
40 {
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
45 {
46  const fvMesh& mesh = ssf.mesh();
47 
48  const labelUList& owner = mesh.owner();
49  const labelUList& neighbour = mesh.neighbour();
50 
51  const volVectorField& C = mesh.C();
52  const surfaceVectorField& Cf = mesh.Cf();
53  const surfaceVectorField& Sf = mesh.Sf();
54  const surfaceScalarField& magSf = mesh.magSf();
55 
56  tmp<volScalarField> treconField
57  (
59  (
60  "reconstruct("+ssf.name()+')',
61  mesh,
63  (
64  ssf.dimensions()/dimArea,
65  scalar(0)
66  ),
67  extrapolatedCalculatedFvPatchScalarField::typeName
68  )
69  );
70  scalarField& rf = treconField.ref();
71 
72  forAll(owner, facei)
73  {
74  label own = owner[facei];
75  label nei = neighbour[facei];
76 
77  rf[own] += (Sf[facei] & (Cf[facei] - C[own]))*ssf[facei]/magSf[facei];
78  rf[nei] -= (Sf[facei] & (Cf[facei] - C[nei]))*ssf[facei]/magSf[facei];
79  }
80 
82 
83  forAll(bsf, patchi)
84  {
85  const fvsPatchScalarField& psf = bsf[patchi];
86 
87  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
88  const vectorField& pCf = Cf.boundaryField()[patchi];
89  const vectorField& pSf = Sf.boundaryField()[patchi];
90  const scalarField& pMagSf = magSf.boundaryField()[patchi];
91 
92  forAll(pOwner, pFacei)
93  {
94  label own = pOwner[pFacei];
95  rf[own] +=
96  (pSf[pFacei] & (pCf[pFacei] - C[own]))
97  *psf[pFacei]/pMagSf[pFacei];
98  }
99  }
100 
101  rf /= mesh.V();
102 
103  treconField.ref().correctBoundaryConditions();
104 
105  return treconField;
106 }
107 
108 
110 {
112  (
113  fvc::reconstructMag(tssf())
114  );
115  tssf.clear();
116  return tvf;
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 
122 } // End namespace fvc
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 } // End namespace Foam
127 
128 // ************************************************************************* //
Foam::surfaceFields.
Graphite solid properties.
Definition: C.H:48
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const surfaceVectorField & Sf() const
Return cell face area vectors.
const dimensionSet dimArea
const word & name() const
Return name.
Definition: IOobject.H:315
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
const surfaceVectorField & Cf() const
Return face centres.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Reconstruct volField from a face flux field.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Generic GeometricField class.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
fvMesh & mesh
tmp< volScalarField > reconstructMag(const surfaceScalarField &)
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:423
const dimensionSet & dimensions() const
Return dimensions.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
const Mesh & mesh() const
Return mesh.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:417
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
const volVectorField & C() const
Return cell centres.
A class for managing temporary objects.
Definition: PtrList.H:53
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Namespace for OpenFOAM.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:800