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 | Copyright (C) 2013 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  (
58  new volScalarField
59  (
60  IOobject
61  (
62  "reconstruct("+ssf.name()+')',
63  ssf.instance(),
64  mesh,
67  ),
68  mesh,
70  (
71  "0",
72  ssf.dimensions()/dimArea,
73  scalar(0)
74  ),
75  zeroGradientFvPatchScalarField::typeName
76  )
77  );
78 
79  scalarField& rf = treconField();
80 
81  forAll(owner, facei)
82  {
83  label own = owner[facei];
84  label nei = neighbour[facei];
85 
86  rf[own] += (Sf[facei] & (Cf[facei] - C[own]))*ssf[facei]/magSf[facei];
87  rf[nei] -= (Sf[facei] & (Cf[facei] - C[nei]))*ssf[facei]/magSf[facei];
88  }
89 
90  const surfaceScalarField::GeometricBoundaryField& bsf = ssf.boundaryField();
91 
92  forAll(bsf, patchi)
93  {
94  const fvsPatchScalarField& psf = bsf[patchi];
95 
96  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
97  const vectorField& pCf = Cf.boundaryField()[patchi];
98  const vectorField& pSf = Sf.boundaryField()[patchi];
99  const scalarField& pMagSf = magSf.boundaryField()[patchi];
100 
101  forAll(pOwner, pFacei)
102  {
103  label own = pOwner[pFacei];
104  rf[own] +=
105  (pSf[pFacei] & (pCf[pFacei] - C[own]))
106  *psf[pFacei]/pMagSf[pFacei];
107  }
108  }
109 
110  rf /= mesh.V();
111 
112  treconField().correctBoundaryConditions();
113 
114  return treconField;
115 }
116 
117 
119 {
121  (
122  fvc::reconstructMag(tssf())
123  );
124  tssf.clear();
125  return tvf;
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130 
131 } // End namespace fvc
132 
133 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134 
135 } // End namespace Foam
136 
137 // ************************************************************************* //
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::surfaceFields.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
const surfaceVectorField & Sf() const
Return cell face area vectors.
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:174
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
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:59
const fileName & instance() const
Definition: IOobject.H:337
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:282
Graphite solid properties.
Definition: C.H:57
const Mesh & mesh() const
Return mesh.
dynamicFvMesh & mesh
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
Namespace for OpenFOAM.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const volVectorField & C() const
Return cell centres as volVectorField.
Reconstruct volField from a face flux field.
#define forAll(list, i)
Definition: UList.H:421
label patchi
const dimensionSet & dimensions() const
Return dimensions.
const word & name() const
Return name.
Definition: IOobject.H:260
tmp< volScalarField > reconstructMag(const surfaceScalarField &)
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 labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:288
A class for managing temporary objects.
Definition: PtrList.H:118
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:552