linearUpwind.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) 2011-2012 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 "linearUpwind.H"
27 #include "fvMesh.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 template<class Type>
34 (
36 ) const
37 {
38  const fvMesh& mesh = this->mesh();
39 
41  (
43  (
44  IOobject
45  (
46  "linearUpwind::correction(" + vf.name() + ')',
47  mesh.time().timeName(),
48  mesh,
49  IOobject::NO_READ,
50  IOobject::NO_WRITE,
51  false
52  ),
53  mesh,
54  dimensioned<Type>(vf.name(), vf.dimensions(), pTraits<Type>::zero)
55  )
56  );
57 
59 
60  const surfaceScalarField& faceFlux = this->faceFlux_;
61 
62  const labelList& owner = mesh.owner();
63  const labelList& neighbour = mesh.neighbour();
64 
65  const volVectorField& C = mesh.C();
66  const surfaceVectorField& Cf = mesh.Cf();
67 
68  tmp
69  <
71  <
74  volMesh
75  >
76  > tgradVf = gradScheme_().grad(vf, gradSchemeName_);
77 
78  const GeometricField
79  <
82  volMesh
83  >& gradVf = tgradVf();
84 
85  forAll(faceFlux, facei)
86  {
87  label celli = (faceFlux[facei] > 0) ? owner[facei] : neighbour[facei];
88  sfCorr[facei] = (Cf[facei] - C[celli]) & gradVf[celli];
89  }
90 
91 
93  GeometricBoundaryField& bSfCorr = sfCorr.boundaryField();
94 
95  forAll(bSfCorr, patchi)
96  {
97  fvsPatchField<Type>& pSfCorr = bSfCorr[patchi];
98 
99  if (pSfCorr.coupled())
100  {
101  const labelUList& pOwner =
102  mesh.boundary()[patchi].faceCells();
103 
104  const vectorField& pCf = Cf.boundaryField()[patchi];
105 
106  const scalarField& pFaceFlux = faceFlux.boundaryField()[patchi];
107 
109  (
110  gradVf.boundaryField()[patchi].patchNeighbourField()
111  );
112 
113  // Build the d-vectors
114  vectorField pd(Cf.boundaryField()[patchi].patch().delta());
115 
116  forAll(pOwner, facei)
117  {
118  label own = pOwner[facei];
119 
120  if (pFaceFlux[facei] > 0)
121  {
122  pSfCorr[facei] = (pCf[facei] - C[own]) & gradVf[own];
123  }
124  else
125  {
126  pSfCorr[facei] =
127  (pCf[facei] - pd[facei] - C[own]) & pGradVfNei[facei];
128  }
129  }
130  }
131  }
132 
133  return tsfCorr;
134 }
135 
136 
137 namespace Foam
138 {
139  //makelimitedSurfaceInterpolationScheme(linearUpwind)
142 }
143 
144 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::linearUpwind
linearUpwind interpolation scheme class derived from upwind and returns upwind weighting factors and ...
Definition: linearUpwind.H:52
Foam::volMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:47
Foam::outerProduct::type
typeOfRank< typename pTraits< arg1 >::cmptType, int(pTraits< arg1 >::rank)+int(pTraits< arg2 >::rank) >::type type
Definition: products.H:72
Foam::fvsPatchField< Type >
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
linearUpwind.H
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::linearUpwind::correction
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the face-interpolate.
Definition: linearUpwind.C:34
makelimitedSurfaceInterpolationTypeScheme
#define makelimitedSurfaceInterpolationTypeScheme(SS, Type)
Definition: limitedSurfaceInterpolationScheme.H:202
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned< Type >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::UList
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
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::C
Graphite solid properties.
Definition: C.H:57
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvsPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvsPatchField.H:308