extendedCellToCellStencilTemplates.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 -------------------------------------------------------------------------------
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 
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
33 template<class Type, class WeightType>
35 <
37  <
41  >
43 (
44  const mapDistribute& map,
45  const labelListList& stencil,
46  const GeometricField<Type, fvPatchField, volMesh>& fld,
47  const List<List<WeightType>>& stencilWeights
48 )
49 {
50  typedef typename outerProduct<WeightType, Type>::type WeightedType;
51  typedef GeometricField<WeightedType, fvPatchField, volMesh>
52  WeightedFieldType;
53 
54  const fvMesh& mesh = fld.mesh();
55 
56  // Collect internal and boundary values
57  List<List<Type>> stencilFld;
58  extendedCellToFaceStencil::collectData(map, stencil, fld, stencilFld);
59 
60  tmp<WeightedFieldType> twf
61  (
62  new WeightedFieldType
63  (
64  IOobject
65  (
66  fld.name(),
67  mesh.time().timeName(),
68  mesh
69  ),
70  mesh,
71  dimensioned<WeightedType>(fld.dimensions(), Zero)
72  )
73  );
74  WeightedFieldType& wf = twf();
75 
76  forAll(wf, celli)
77  {
78  const List<Type>& stField = stencilFld[celli];
79  const List<WeightType>& stWeight = stencilWeights[celli];
80 
81  forAll(stField, i)
82  {
83  wf[celli] += stWeight[i]*stField[i];
84  }
85  }
86 
87  // Boundaries values?
88 
89  return twf;
90 }
91 
92 
93 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
extendedCellToFaceStencil.H
Foam::extendedCellToFaceStencil::collectData
static void collectData(const mapDistribute &map, const labelListList &stencil, const GeometricField< T, fvPatchField, volMesh > &fld, List< List< T >> &stencilFld)
Use map to get the data into stencil order.
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::volMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:51
Foam::outerProduct::type
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:114
Foam::extendedCellToCellStencil::weightedSum
static tmp< GeometricField< typename outerProduct< WeightType, Type >::type, fvPatchField, volMesh > > weightedSum(const mapDistribute &map, const labelListList &stencil, const GeometricField< Type, fvPatchField, volMesh > &fld, const List< List< WeightType >> &stencilWeights)
Sum surface field contributions to create cell values.
extendedCellToCellStencil.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
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::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53