volPointInterpolationAdjoint.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  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "fvMesh.H"
31 #include "volFields.H"
32 #include "pointFields.H"
33 #include "pointConstraints.H"
34 #include "surfaceFields.H"
35 #include "processorPointPatch.H"
36 #include "symmetryPolyPatch.H"
37 #include "symmetryPlanePolyPatch.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
50 {
51  if (debug)
52  {
53  Pout<< "volPointInterpolationAdjoint::calcBoundaryAddressing() : "
54  << "constructing boundary addressing"
55  << endl;
56  }
57 
58  boundaryPtr_.reset
59  (
60  new primitivePatch
61  (
62  SubList<face>
63  (
64  mesh().faces(),
65  mesh().nBoundaryFaces(),
66  mesh().nInternalFaces()
67  ),
68  mesh().points()
69  )
70  );
72 
74  boundaryIsPatchFace_ = false;
75 
76  // Store per mesh point whether it is on any 'real' patch. Currently
77  // boolList just so we can use syncUntransformedData (does not take
78  // bitSet. Tbd)
79  boolList isPatchPoint(mesh().nPoints(), false);
80  boolList isSymmetryPoint(mesh().nPoints(), false);
81 
82  const polyBoundaryMesh& pbm = mesh().boundaryMesh();
83 
84  // Get precalculated volField only so we can use coupled() tests for
85  // cyclicAMI
86  const surfaceScalarField& magSf = mesh().magSf();
87 
88  forAll(pbm, patchi)
89  {
90  const polyPatch& pp = pbm[patchi];
91 
92  if
93  (
94  !isA<emptyPolyPatch>(pp)
95  && !magSf.boundaryField()[patchi].coupled()
96  && !isA<symmetryPolyPatch>(pp)
97  && !isA<symmetryPlanePolyPatch>(pp)
98  )
99  {
100  label bFacei = pp.start()-mesh().nInternalFaces();
101 
102  forAll(pp, i)
103  {
104  boundaryIsPatchFace_[bFacei] = true;
105 
106  const face& f = boundary[bFacei++];
107 
108  forAll(f, fp)
109  {
110  isPatchPoint[f[fp]] = true;
111  }
112  }
113  }
114  else if (isA<symmetryPolyPatch>(pp) || isA<symmetryPlanePolyPatch>(pp))
115  {
116  const labelList& meshPoints = pp.meshPoints();
117  for (const label pointI : meshPoints)
118  {
119  isSymmetryPoint[pointI] = true;
120  }
121  }
122  }
123 
124  // Make sure point status is synchronised so even processor that holds
125  // no face of a certain patch still can have boundary points marked.
127  (
128  mesh(),
129  isPatchPoint,
130  orEqOp<bool>()
131  );
132 
133  // Convert to bitSet
135  isPatchPoint_.assign(isPatchPoint);
136 
138  isSymmetryPoint_.assign(isSymmetryPoint);
139 
140  if (debug)
141  {
142  label nPatchFace = 0;
144  {
145  if (boundaryIsPatchFace_[i])
146  {
147  nPatchFace++;
148  }
149  }
150  label nPatchPoint = 0;
152  {
153  if (isPatchPoint_[i])
154  {
155  nPatchPoint++;
156  }
157  }
158  Pout<< "boundary:" << nl
159  << " faces :" << boundary.size() << nl
160  << " of which on proper patch:" << nPatchFace << nl
161  << " points:" << boundary.nPoints() << nl
162  << " of which on proper patch:" << nPatchPoint << endl;
163  }
164 }
165 
166 
168 (
169  scalarField& sumWeights
170 )
171 {
172  if (debug)
173  {
174  Pout<< "volPointInterpolationAdjoint::makeBoundaryWeights() : "
175  << "constructing weighting factors for boundary points." << endl;
176  }
177 
178  const pointField& points = mesh().points();
179  const pointField& faceCentres = mesh().faceCentres();
180 
181  const primitivePatch& boundary = boundaryPtr_();
182 
183  boundaryPointWeights_.clear();
184  boundaryPointWeights_.setSize(boundary.meshPoints().size());
185 
186  forAll(boundary.meshPoints(), i)
187  {
188  label pointi = boundary.meshPoints()[i];
189 
190  if (isPatchPoint_[pointi])
191  {
192  const labelList& pFaces = boundary.pointFaces()[i];
193 
194  scalarList& pw = boundaryPointWeights_[i];
195  pw.setSize(pFaces.size());
196 
197  sumWeights[pointi] = 0.0;
198 
199  forAll(pFaces, i)
200  {
201  if (boundaryIsPatchFace_[pFaces[i]])
202  {
203  label facei = mesh().nInternalFaces() + pFaces[i];
204 
205  pw[i] = 1.0/mag(points[pointi] - faceCentres[facei]);
206  sumWeights[pointi] += pw[i];
207  }
208  else
209  {
210  pw[i] = 0.0;
211  }
212  }
213  }
214  }
215 }
216 
217 
219 {
220  if (debug)
221  {
222  Pout<< "volPointInterpolationAdjoint::makeWeights() : "
223  << "constructing weighting factors"
224  << endl;
225  }
226 
227  const pointMesh& pMesh = pointMesh::New(mesh());
228 
229  // Update addressing over all boundary faces
230  calcBoundaryAddressing();
231 
232 
233  // Running sum of weights
234  tmp<pointScalarField> tsumWeights
235  (
236  new pointScalarField
237  (
238  IOobject
239  (
240  "volPointSumWeights",
242  mesh()
243  ),
244  pMesh,
246  )
247  );
248  pointScalarField& sumWeights = tsumWeights.ref();
249 
250 
251  // Create boundary weights; sumWeights
252  makeBoundaryWeights(sumWeights);
253 
254 
255  const primitivePatch& boundary = boundaryPtr_();
256  const labelList& mp = boundary.meshPoints();
257 
258 
259  // Sum collocated contributions
261  (
262  mesh(),
263  sumWeights,
264  plusEqOp<scalar>()
265  );
266 
267 
268  // Push master data to slaves. It is possible (not sure how often) for
269  // a coupled point to have its master on a different patch so
270  // to make sure just push master data to slaves. Reuse the syncPointData
271  // structure.
272  pushUntransformedData(sumWeights);
273 
274  // Normalise boundary weights
275  forAll(mp, i)
276  {
277  const label pointi = mp[i];
278 
279  scalarList& pw = boundaryPointWeights_[i];
280  // Note:pw only sized for isPatchPoint
281  forAll(pw, i)
282  {
283  pw[i] /= sumWeights[pointi];
284  }
285  }
286 
287  if (debug)
288  {
289  Pout<< "volPointInterpolationAdjoint::makeWeights() : "
290  << "finished constructing weighting factors"
291  << endl;
292  }
293 }
294 
295 
296 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
297 
299 :
300  MeshObject<fvMesh, Foam::UpdateableMeshObject, volPointInterpolationAdjoint>
301  (
302  vm
303  )
304 {
305  makeWeights();
306 }
307 
308 
309 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * //
312 {}
313 
314 
315 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
318 {
319  makeWeights();
320 }
321 
322 
324 {
325  makeWeights();
326 
327  return true;
328 }
329 
330 
332 (
333  const vectorField& pf,
334  vectorField& vf,
335  const labelHashSet& patchIDs
336 ) const
337 {
338  // Get face data in flat list
339  const fvMesh& Mesh = mesh();
340 
341  vectorField boundaryVals(Mesh.nBoundaryFaces(), Zero);
342 
343  // Do points on 'normal' patches from the surrounding patch faces
344  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
345 
346  const primitivePatch& boundary = boundaryPtr_();
347  const labelList& mp = boundary.meshPoints();
348 
349  forAll(mp, i)
350  {
351  label pointi = mp[i];
352 
353  if (isPatchPoint_[pointi])
354  {
355  const labelList& pFaces = boundary.pointFaces()[i];
356  const scalarList& pWeights = boundaryPointWeights_[i];
357  const vector& val = pf[pointi];
358  // In symmetry planes, face-to-point weights should, in general,
359  // have half the weight of what they actually do in
360  // volPointInterpolation since, in a complete case, a face laying
361  // on the opposite side of the symmetry plane would also contribute
362  // to a point laying on the symmetry plane.
363  // For face-to-point interpolation this is not a problem, but for
364  // the adjoint point-to-face interpolation, the correct value of
365  // the weight should be taken into consideration
366  scalar mod(isSymmetryPoint_[pointi] ? 0.5 : 1);
367 
368  forAll(pFaces, j)
369  {
370  if (boundaryIsPatchFace_[pFaces[j]])
371  {
372  boundaryVals[pFaces[j]] += mod*pWeights[j]*val;
373  }
374  }
375  }
376  }
377 
378  // Transfer values to face-based sensitivity field
379  label nPassedFaces(0);
380  for (const label patchi : patchIDs)
381  {
382  const fvPatch& patch = Mesh.boundary()[patchi];
383  label bFacei = patch.start() - Mesh.nInternalFaces();
384  SubList<vector> patchFaceSens(vf, patch.size(), nPassedFaces);
385  patchFaceSens = SubList<vector>(boundaryVals, patch.size(), bFacei);
386  nPassedFaces += patch.size();
387  }
388 }
389 
390 
391 // ************************************************************************* //
Foam::surfaceFields.
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
faceListList boundary
const labelList patchIDs(pbm.indices(polyPatchNames, true))
const polyBoundaryMesh & pbm
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
bitSet isPatchPoint_
Per mesh(!) point whether is on non-coupled, non-empty patch (on.
void calcBoundaryAddressing()
Construct addressing over all boundary faces.
void updateMesh(const mapPolyMesh &)
Update mesh topology using the morph engine.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
void setSize(const label n, unsigned int val=0u)
Alias for resize()
Definition: PackedList.H:785
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const dimensionSet dimless
Dimensionless.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
bitSet isSymmetryPoint_
Per mesh(!) point whether is on symmetry plane.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
autoPtr< primitivePatch > boundaryPtr_
Boundary addressing.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Interpolate from cell centres to points (vertices) using inverse distance weighting.
A list of faces which address into the list of points.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
void makeWeights()
Construct all point weighting factors.
void setSize(const label n)
Alias for resize()
Definition: List.H:316
bool movePoints()
Correct weighting factors for moving mesh.
dynamicFvMesh & mesh
const pointField & points
const fvMesh & mesh() const noexcept
Reference to the mesh.
Definition: MeshObject.H:157
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static void syncUntransformedData(const polyMesh &mesh, List< Type > &pointData, const CombineOp &cop)
Helper: sync data on collocated points only.
bitSet boundaryIsPatchFace_
Per boundary face whether is on non-coupled, non-empty patch.
label nInternalFaces() const noexcept
Number of internal faces.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field...
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:266
const vectorField & faceCentres() const
void makeBoundaryWeights(scalarField &sumWeights)
Make weights for points on uncoupled patches.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
volPointInterpolationAdjoint(const volPointInterpolationAdjoint &)=delete
No copy construct.
const std::string patch
OpenFOAM patch number as a std::string.
void assign(const UList< bool > &bools)
Copy assign all entries from a list of bools.
Definition: bitSet.C:267
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
List< label > labelList
A List of labels.
Definition: List.H:62
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
void interpolateSensitivitiesField(const GeometricField< Type, pointPatchField, pointMesh > &pf, typename GeometricField< Type, fvPatchField, volMesh >::Boundary &vf, const labelHashSet &patchIDs) const
Interpolate sensitivties from points to faces.
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
List< bool > boolList
A List of bools.
Definition: List.H:60
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const dimensionedScalar mp
Proton mass.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127