extendedUpwindCellToFaceStencilTemplates.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) 2011-2016 OpenFOAM Foundation
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 
29 
30 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
31 
32 template<class Type>
35 (
36  const surfaceScalarField& phi,
38  const List<List<scalar>>& ownWeights,
39  const List<List<scalar>>& neiWeights
40 ) const
41 {
42  const fvMesh& mesh = fld.mesh();
43 
44  // Collect internal and boundary values
45  List<List<Type>> ownFld;
46  collectData(ownMap(), ownStencil(), fld, ownFld);
47  List<List<Type>> neiFld;
48  collectData(neiMap(), neiStencil(), fld, neiFld);
49 
51  (
53  (
54  IOobject
55  (
56  fld.name(),
57  mesh.time().timeName(),
58  mesh,
59  IOobject::NO_READ,
60  IOobject::NO_WRITE,
61  false
62  ),
63  mesh,
64  dimensioned<Type>(fld.dimensions(), Zero)
65  )
66  );
68 
69  // Internal faces
70  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
71  {
72  if (phi[facei] > 0)
73  {
74  // Flux out of owner. Use upwind (= owner side) stencil.
75  const List<Type>& stField = ownFld[facei];
76  const List<scalar>& stWeight = ownWeights[facei];
77 
78  forAll(stField, i)
79  {
80  sf[facei] += stField[i]*stWeight[i];
81  }
82  }
83  else
84  {
85  const List<Type>& stField = neiFld[facei];
86  const List<scalar>& stWeight = neiWeights[facei];
87 
88  forAll(stField, i)
89  {
90  sf[facei] += stField[i]*stWeight[i];
91  }
92  }
93  }
94 
95  // Boundaries. Either constrained or calculated so assign value
96  // directly (instead of nicely using operator==)
98  Boundary& bSfCorr = sf.boundaryFieldRef();
99 
100  forAll(bSfCorr, patchi)
101  {
102  fvsPatchField<Type>& pSfCorr = bSfCorr[patchi];
103 
104  if (pSfCorr.coupled())
105  {
106  label facei = pSfCorr.patch().start();
107 
108  forAll(pSfCorr, i)
109  {
110  if (phi.boundaryField()[patchi][i] > 0)
111  {
112  // Flux out of owner. Use upwind (= owner side) stencil.
113  const List<Type>& stField = ownFld[facei];
114  const List<scalar>& stWeight = ownWeights[facei];
115 
116  forAll(stField, j)
117  {
118  pSfCorr[i] += stField[j]*stWeight[j];
119  }
120  }
121  else
122  {
123  const List<Type>& stField = neiFld[facei];
124  const List<scalar>& stWeight = neiWeights[facei];
125 
126  forAll(stField, j)
127  {
128  pSfCorr[i] += stField[j]*stWeight[j];
129  }
130  }
131  facei++;
132  }
133  }
134  }
135 
136  return tsfCorr;
137 }
138 
139 
140 // ************************************************************************* //
extendedCellToFaceStencil.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::extendedUpwindCellToFaceStencil::weightedSum
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > weightedSum(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &fld, const List< List< scalar >> &ownWeights, const List< List< scalar >> &neiWeights) const
Sum vol field contributions to create face values.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::fvsPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvsPatchField.H:281
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:749
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::fvsPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvsPatchField.H:310