extendedFaceToCellStencilTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 distributionMap& map,
34  const labelListList& stencil,
35  const SurfaceField<Type>& fld,
36  List<List<Type>>& stencilFld
37 )
38 {
39  // 1. Construct face data in compact addressing
40  List<Type> flatFld(map.constructSize(), Zero);
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 fvsPatchField<Type>& pfld = fld.boundaryField()[patchi];
51 
52  label nCompact = pfld.patch().start();
53 
54  forAll(pfld, i)
55  {
56  flatFld[nCompact++] = pfld[i];
57  }
58  }
59 
60  // Do all swapping
61  map.distribute(flatFld);
62 
63  // 2. Pull to stencil
64  stencilFld.setSize(stencil.size());
65 
66  forAll(stencil, facei)
67  {
68  const labelList& compactCells = stencil[facei];
69 
70  stencilFld[facei].setSize(compactCells.size());
71 
72  forAll(compactCells, i)
73  {
74  stencilFld[facei][i] = flatFld[compactCells[i]];
75  }
76  }
77 }
78 
79 
80 template<class Type>
83 (
84  const distributionMap& map,
85  const labelListList& stencil,
86  const SurfaceField<Type>& fld,
87  const List<List<scalar>>& stencilWeights
88 )
89 {
90  const fvMesh& mesh = fld.mesh();
91 
92  // Collect internal and boundary values
93  List<List<Type>> stencilFld;
94  collectData(map, stencil, fld, stencilFld);
95 
96  tmp<VolField<Type>> tsfCorr
97  (
99  (
100  fld.name(),
101  mesh,
103  (
104  fld.name(),
105  fld.dimensions(),
106  Zero
107  )
108  )
109  );
110  VolField<Type>& sf = tsfCorr.ref();
111 
112  // cells
113  forAll(sf, celli)
114  {
115  const List<Type>& stField = stencilFld[celli];
116  const List<scalar>& stWeight = stencilWeights[celli];
117 
118  forAll(stField, i)
119  {
120  sf[celli] += stField[i]*stWeight[i];
121  }
122  }
123 
124  // Boundaries values?
125 
126  return tsfCorr;
127 }
128 
129 
130 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
Generic dimensioned Type class.
label constructSize() const
Constructed data size.
Class containing processor-to-processor mapping information.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
static void collectData(const distributionMap &map, const labelListList &stencil, const SurfaceField< Type > &fld, List< List< Type >> &stencilFld)
Use map to get the data into stencil order.
static tmp< VolField< Type > > weightedSum(const distributionMap &map, const labelListList &stencil, const SurfaceField< Type > &fld, const List< List< scalar >> &stencilWeights)
Sum surface field contributions to create cell values.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:151
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
const fvPatch & patch() const
Return patch.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
label patchi
volScalarField sf(fieldObject, mesh)
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
static const zero Zero
Definition: zero.H:97
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