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 | 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>
34 (
35  const mapDistribute& map,
36  const labelListList& stencil,
38  List<List<Type>>& stencilFld
39 )
40 {
41  // 1. Construct cell data in compact addressing
42  List<Type> flatFld(map.constructSize(), Zero);
43 
44  // Insert my internal values
45  forAll(fld, celli)
46  {
47  flatFld[celli] = fld[celli];
48  }
49  // Insert my boundary values
50  forAll(fld.boundaryField(), patchi)
51  {
52  const fvPatchField<Type>& pfld = fld.boundaryField()[patchi];
53 
54  label nCompact =
55  pfld.patch().start()
56  -fld.mesh().nInternalFaces()
57  +fld.mesh().nCells();
58 
59  forAll(pfld, i)
60  {
61  flatFld[nCompact++] = pfld[i];
62  }
63  }
64 
65  // Do all swapping
66  map.distribute(flatFld);
67 
68  // 2. Pull to stencil
69  stencilFld.setSize(stencil.size());
70 
71  forAll(stencil, facei)
72  {
73  const labelList& compactCells = stencil[facei];
74 
75  stencilFld[facei].setSize(compactCells.size());
76 
77  forAll(compactCells, i)
78  {
79  stencilFld[facei][i] = flatFld[compactCells[i]];
80  }
81  }
82 }
83 
84 
85 template<class Type>
88 (
89  const mapDistribute& map,
90  const labelListList& stencil,
91  const GeometricField<Type, fvPatchField, volMesh>& fld,
92  const List<List<scalar>>& stencilWeights
93 )
94 {
95  const fvMesh& mesh = fld.mesh();
96 
97  // Collect internal and boundary values
98  List<List<Type>> stencilFld;
99  collectData(map, stencil, fld, stencilFld);
100 
101  tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tsfCorr
102  (
103  new GeometricField<Type, fvsPatchField, surfaceMesh>
104  (
105  IOobject
106  (
107  fld.name(),
108  mesh.time().timeName(),
109  mesh,
113  ),
114  mesh,
115  dimensioned<Type>(fld.dimensions(), Zero)
116  )
117  );
118  GeometricField<Type, fvsPatchField, surfaceMesh>& sf = tsfCorr.ref();
119 
120  // Internal faces
121  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
122  {
123  const List<Type>& stField = stencilFld[facei];
124  const List<scalar>& stWeight = stencilWeights[facei];
125 
126  forAll(stField, i)
127  {
128  sf[facei] += stField[i]*stWeight[i];
129  }
130  }
131 
132  // Boundaries. Either constrained or calculated so assign value
133  // directly (instead of nicely using operator==)
134  typename GeometricField<Type, fvsPatchField, surfaceMesh>::
135  Boundary& bSfCorr = sf.boundaryFieldRef();
136 
137  forAll(bSfCorr, patchi)
138  {
139  fvsPatchField<Type>& pSfCorr = bSfCorr[patchi];
140 
141  if (pSfCorr.coupled())
142  {
143  label facei = pSfCorr.patch().start();
144 
145  forAll(pSfCorr, i)
146  {
147  const List<Type>& stField = stencilFld[facei];
148  const List<scalar>& stWeight = stencilWeights[facei];
149 
150  forAll(stField, j)
151  {
152  pSfCorr[i] += stField[j]*stWeight[j];
153  }
154 
155  facei++;
156  }
157  }
158  }
159 
160  return tsfCorr;
161 }
162 
163 
164 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute List data using default commsType, default flip/negate operator.
label constructSize() const noexcept
Constructed data size.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
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.
label nInternalFaces() const noexcept
Number of internal faces.
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.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
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))
Class containing processor-to-processor mapping information.
Nothing to be read.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Do not request registration (bool: false)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127