Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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"
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
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 const polyBoundaryMesh& pbm = mesh().boundaryMesh();
59
60 boundaryPtr_.reset(new primitivePatch(pbm.faces(), mesh().points()));
61 const auto& boundary = *boundaryPtr_;
62
65
66 // Store per mesh point whether it is on any 'real' patch. Currently
67 // boolList just so we can use syncUntransformedData (does not take
68 // bitSet. Tbd)
69 boolList isPatchPoint(mesh().nPoints(), false);
70 boolList isSymmetryPoint(mesh().nPoints(), false);
71
72 // Get precalculated volField only so we can use coupled() tests for
73 // cyclicAMI
74 const surfaceScalarField& magSf = mesh().magSf();
75
76 forAll(pbm, patchi)
77 {
78 const polyPatch& pp = pbm[patchi];
79
80 if
81 (
83 && !magSf.boundaryField()[patchi].coupled()
86 )
87 {
88 boundaryIsPatchFace_.set(labelRange(pp.offset(), pp.size()));
89
90 for (const auto& f : pp.faces())
91 {
92 UIndirectList<bool>(isPatchPoint, f) = true;
93 }
94 }
96 {
97 UIndirectList<bool>(isSymmetryPoint, pp.meshPoints()) = true;
98 }
99 }
100
101 // Make sure point status is synchronised so even processor that holds
102 // no face of a certain patch still can have boundary points marked.
104 (
105 mesh(),
106 isPatchPoint,
107 orEqOp<bool>()
108 );
109
110 // Convert to bitSet
111 isPatchPoint_.reset();
112 isPatchPoint_.resize(mesh().nPoints());
113 isPatchPoint_.assign(isPatchPoint);
114
115 isSymmetryPoint_.reset();
116 isSymmetryPoint_.resize(mesh().nPoints());
117 isSymmetryPoint_.assign(isSymmetryPoint);
118
119 if (debug)
120 {
121 label nPatchFace = boundaryIsPatchFace_.count();
122 label nPatchPoint = isPatchPoint_.count();
123
124 Pout<< "boundary:" << nl
125 << " faces :" << boundary.size() << nl
126 << " of which on proper patch:" << nPatchFace << nl
127 << " points:" << boundary.nPoints() << nl
128 << " of which on proper patch:" << nPatchPoint << endl;
129 }
130}
131
132
134(
135 scalarField& sumWeights
136)
137{
138 if (debug)
139 {
140 Pout<< "volPointInterpolationAdjoint::makeBoundaryWeights() : "
141 << "constructing weighting factors for boundary points." << endl;
142 }
143
144 const pointField& points = mesh().points();
145 const pointField& faceCentres = mesh().faceCentres();
146
147 const primitivePatch& boundary = boundaryPtr_();
148
149 boundaryPointWeights_.clear();
150 boundaryPointWeights_.setSize(boundary.meshPoints().size());
151
152 forAll(boundary.meshPoints(), i)
153 {
154 label pointi = boundary.meshPoints()[i];
155
156 if (isPatchPoint_[pointi])
157 {
158 const labelList& pFaces = boundary.pointFaces()[i];
159
160 scalarList& pw = boundaryPointWeights_[i];
161 pw.setSize(pFaces.size());
162
163 sumWeights[pointi] = 0.0;
164
165 forAll(pFaces, i)
166 {
167 if (boundaryIsPatchFace_[pFaces[i]])
168 {
169 label facei = mesh().nInternalFaces() + pFaces[i];
170
171 pw[i] = 1.0/mag(points[pointi] - faceCentres[facei]);
172 sumWeights[pointi] += pw[i];
173 }
174 else
175 {
176 pw[i] = 0.0;
178 }
179 }
180 }
181}
182
183
185{
186 if (debug)
187 {
188 Pout<< "volPointInterpolationAdjoint::makeWeights() : "
189 << "constructing weighting factors"
190 << endl;
191 }
192
193 const pointMesh& pMesh = pointMesh::New(mesh());
194
195 // Update addressing over all boundary faces
196 calcBoundaryAddressing();
197
198
199 // Running sum of weights
200 tmp<pointScalarField> tsumWeights
201 (
203 (
205 (
206 "volPointSumWeights",
208 mesh()
209 ),
210 pMesh,
212 )
213 );
214 pointScalarField& sumWeights = tsumWeights.ref();
215
216
217 // Create boundary weights; sumWeights
218 makeBoundaryWeights(sumWeights);
219
220
221 const primitivePatch& boundary = boundaryPtr_();
222 const labelList& mp = boundary.meshPoints();
223
224
225 // Sum collocated contributions
227 (
228 mesh(),
229 sumWeights,
231 );
232
233
234 // Push master data to slaves. It is possible (not sure how often) for
235 // a coupled point to have its master on a different patch so
236 // to make sure just push master data to slaves. Reuse the syncPointData
237 // structure.
238 pushUntransformedData(sumWeights);
239
240 // Normalise boundary weights
241 forAll(mp, i)
242 {
243 const label pointi = mp[i];
244
245 scalarList& pw = boundaryPointWeights_[i];
246 // Note:pw only sized for isPatchPoint
247 forAll(pw, i)
248 {
249 pw[i] /= sumWeights[pointi];
250 }
251 }
252
253 if (debug)
254 {
255 Pout<< "volPointInterpolationAdjoint::makeWeights() : "
256 << "finished constructing weighting factors"
258 }
259}
260
261
262// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
263
265:
266 MeshObject_type(vm)
268 makeWeights();
269}
270
271
272// * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * //
277
278// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
285
288 makeWeights();
289
290 return true;
291}
292
293
295(
296 const vectorField& pf,
297 vectorField& vf,
299) const
300{
301 // Get face data in flat list
302 const fvMesh& Mesh = mesh();
303
304 vectorField boundaryVals(Mesh.nBoundaryFaces(), Zero);
305
306 // Do points on 'normal' patches from the surrounding patch faces
307 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
308
309 const primitivePatch& boundary = boundaryPtr_();
310 const labelList& mp = boundary.meshPoints();
311
312 forAll(mp, i)
313 {
314 label pointi = mp[i];
315
316 if (isPatchPoint_[pointi])
317 {
318 const labelList& pFaces = boundary.pointFaces()[i];
319 const scalarList& pWeights = boundaryPointWeights_[i];
320 const vector& val = pf[pointi];
321 // In symmetry planes, face-to-point weights should, in general,
322 // have half the weight of what they actually do in
323 // volPointInterpolation since, in a complete case, a face laying
324 // on the opposite side of the symmetry plane would also contribute
325 // to a point laying on the symmetry plane.
326 // For face-to-point interpolation this is not a problem, but for
327 // the adjoint point-to-face interpolation, the correct value of
328 // the weight should be taken into consideration
329 scalar mod(isSymmetryPoint_[pointi] ? 0.5 : 1);
330
331 forAll(pFaces, j)
332 {
333 if (boundaryIsPatchFace_[pFaces[j]])
334 {
335 boundaryVals[pFaces[j]] += mod*pWeights[j]*val;
336 }
337 }
338 }
339 }
340
341 // Transfer values to face-based sensitivity field
342 label nPassedFaces(0);
343 for (const label patchi : patchIDs)
344 {
345 const fvPatch& patch = Mesh.boundary()[patchi];
346 label bFacei = patch.start() - Mesh.nInternalFaces();
347 SubList<vector> patchFaceSens(vf, patch.size(), nPassedFaces);
348 patchFaceSens = SubList<vector>(boundaryVals, patch.size(), bFacei);
349 nPassedFaces += patch.size();
350 }
351}
352
353
354// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
const polyBoundaryMesh & pbm
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
void setSize(label n)
Alias for resize().
Definition List.H:536
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
void resize_nocopy(const label numElem)
Currently identical to resize. Subject to future change (Oct-2021).
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A range or interval of labels defined by a start and a size.
Definition labelRange.H:66
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
static void syncUntransformedData(const polyMesh &mesh, List< Type > &pointData, const CombineOp &cop)
Helper: sync data on collocated points only.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces).
const vectorField & faceCentres() const
label nInternalFaces() const noexcept
Number of internal faces.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Interpolate from cell centres to points (vertices) using inverse distance weighting.
bool movePoints()
Correct weighting factors for moving mesh.
void pushUntransformedData(List< Type > &) const
Helper: push master point data to collocated points.
bitSet isPatchPoint_
Per mesh(!) point whether is on non-coupled, non-empty patch (on.
void makeWeights()
Construct all point weighting factors.
scalarListList boundaryPointWeights_
Per boundary point the weights per pointFaces.
volPointInterpolationAdjoint(const volPointInterpolationAdjoint &)=delete
No copy construct.
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.
void calcBoundaryAddressing()
Construct addressing over all boundary faces.
void makeBoundaryWeights(scalarField &sumWeights)
Make weights for points on uncoupled patches.
autoPtr< primitivePatch > boundaryPtr_
Boundary addressing.
void updateMesh(const mapPolyMesh &)
Update mesh topology using the morph engine.
bitSet boundaryIsPatchFace_
Per boundary face whether is on non-coupled, non-empty patch.
bitSet isSymmetryPoint_
Per mesh(!) point whether is on symmetry plane.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
faceListList boundary
dynamicFvMesh & mesh
const pointField & points
label nPoints
const dimensionedScalar mp
Proton mass.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
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];} pointMap[start]=pointMap[end]=Foam::min(start, end);} } 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};constexpr 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< DynamicList< face > > pFaces[nBCs]
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.