extendedCellToFaceStencilTemplates.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-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 
27 
28 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
30 template<class Type>
32 (
33  const mapDistribute& map,
34  const labelListList& stencil,
36  List<List<Type> >& stencilFld
37 )
38 {
39  // 1. Construct cell data in compact addressing
41 
42  // Insert my internal values
43  forAll(fld, cellI)
44  {
45  flatFld[cellI] = fld[cellI];
46  }
47  // Insert my boundary values
48  forAll(fld.boundaryField(), patchI)
49  {
50  const fvPatchField<Type>& pfld = fld.boundaryField()[patchI];
51 
52  label nCompact =
53  pfld.patch().start()
54  -fld.mesh().nInternalFaces()
55  +fld.mesh().nCells();
56 
57  forAll(pfld, i)
58  {
59  flatFld[nCompact++] = pfld[i];
60  }
61  }
62 
63  // Do all swapping
64  map.distribute(flatFld);
65 
66  // 2. Pull to stencil
67  stencilFld.setSize(stencil.size());
68 
69  forAll(stencil, faceI)
70  {
71  const labelList& compactCells = stencil[faceI];
72 
73  stencilFld[faceI].setSize(compactCells.size());
74 
75  forAll(compactCells, i)
76  {
77  stencilFld[faceI][i] = flatFld[compactCells[i]];
78  }
79  }
80 }
81 
82 
83 template<class Type>
86 (
87  const mapDistribute& map,
88  const labelListList& stencil,
90  const List<List<scalar> >& stencilWeights
91 )
92 {
93  const fvMesh& mesh = fld.mesh();
94 
95  // Collect internal and boundary values
96  List<List<Type> > stencilFld;
97  collectData(map, stencil, fld, stencilFld);
98 
100  (
102  (
103  IOobject
104  (
105  fld.name(),
106  mesh.time().timeName(),
107  mesh,
108  IOobject::NO_READ,
109  IOobject::NO_WRITE,
110  false
111  ),
112  mesh,
114  (
115  fld.name(),
116  fld.dimensions(),
118  )
119  )
120  );
122 
123  // Internal faces
124  for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
125  {
126  const List<Type>& stField = stencilFld[faceI];
127  const List<scalar>& stWeight = stencilWeights[faceI];
128 
129  forAll(stField, i)
130  {
131  sf[faceI] += stField[i]*stWeight[i];
132  }
133  }
134 
135  // Boundaries. Either constrained or calculated so assign value
136  // directly (instead of nicely using operator==)
138  GeometricBoundaryField& bSfCorr = sf.boundaryField();
139 
140  forAll(bSfCorr, patchi)
141  {
142  fvsPatchField<Type>& pSfCorr = bSfCorr[patchi];
143 
144  if (pSfCorr.coupled())
145  {
146  label faceI = pSfCorr.patch().start();
147 
148  forAll(pSfCorr, i)
149  {
150  const List<Type>& stField = stencilFld[faceI];
151  const List<scalar>& stWeight = stencilWeights[faceI];
152 
153  forAll(stField, j)
154  {
155  pSfCorr[i] += stField[j]*stWeight[j];
156  }
157 
158  faceI++;
159  }
160  }
161  }
162 
163  return tsfCorr;
164 }
165 
166 
167 // ************************************************************************* //
Foam::fvPatchField< Type >
extendedCellToFaceStencil.H
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::fvsPatchField< Type >
Foam::extendedCellToFaceStencil::weightedSum
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > weightedSum(const mapDistribute &map, const labelListList &stencil, const GeometricField< Type, fvPatchField, volMesh > &fld, const List< List< scalar > > &stencilWeights)
Sum vol field contributions to create face values.
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::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:152
Foam::mapDistribute::distribute
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Definition: mapDistributeTemplates.C:155
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::dimensioned< Type >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
sf
volScalarField sf(fieldObject, mesh)
Foam::mapDistributeBase::constructSize
label constructSize() const
Constructed data size.
Definition: mapDistributeBase.H:244
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::fvsPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvsPatchField.H:278
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
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::fvPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:300
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
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::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