Loading...
Searching...
No Matches
volPointInterpolateAdjoint.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-2017 OpenFOAM Foundation
9 Copyright (C) 2016-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 "volFields.H"
31#include "pointFields.H"
32#include "emptyFvPatch.H"
35#include "symmetryPolyPatch.H"
37
38// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39
40template<class Type>
42(
43 List<Type>& pointData
44) const
45{
46 // Transfer onto coupled patch
47 const globalMeshData& gmd = mesh().globalData();
48 const indirectPrimitivePatch& cpp = gmd.coupledPatch();
49 const labelList& meshPoints = cpp.meshPoints();
50
51 const mapDistribute& slavesMap = gmd.globalCoPointSlavesMap();
52 const labelListList& slaves = gmd.globalCoPointSlaves();
53
54 List<Type> elems(slavesMap.constructSize());
55 forAll(meshPoints, i)
56 {
57 elems[i] = pointData[meshPoints[i]];
58 }
59
60 // Combine master data with slave data
61 forAll(slaves, i)
62 {
63 const labelList& slavePoints = slaves[i];
64
65 // Copy master data to slave slots
66 forAll(slavePoints, j)
67 {
68 elems[slavePoints[j]] = elems[i];
69 }
70 }
71
72 // Push slave-slot data back to slaves
73 slavesMap.reverseDistribute(elems.size(), elems, false);
74
75 // Extract back onto mesh
76 forAll(meshPoints, i)
77 {
78 pointData[meshPoints[i]] = elems[i];
79 }
80}
81
82
83template<class Type>
86(
87 const GeometricField<Type, fvPatchField, volMesh>& vf
88) const
89{
90 const fvMesh& mesh = vf.mesh();
91 const fvBoundaryMesh& bm = mesh.boundary();
92
93 auto tboundaryVals = tmp<Field<Type>>::New(mesh.nBoundaryFaces());
94 auto& boundaryVals = tboundaryVals.ref();
95
96 forAll(vf.boundaryField(), patchi)
97 {
98 label bFacei = bm[patchi].patch().start() - mesh.nInternalFaces();
99
100 if
101 (
102 !isA<emptyFvPatch>(bm[patchi])
103 && !vf.boundaryField()[patchi].coupled()
104 )
105 {
106 SubList<Type>
107 (
108 boundaryVals,
109 vf.boundaryField()[patchi].size(),
110 bFacei
111 ) = vf.boundaryField()[patchi];
112 }
113 else
114 {
115 const polyPatch& pp = bm[patchi].patch();
116
117 forAll(pp, i)
118 {
119 boundaryVals[bFacei++] = Zero;
120 }
121 }
123
124 return tboundaryVals;
125}
126
127
128template<class Type>
130(
131 GeometricField<Type, pointPatchField, pointMesh>& pf
132) const
133{
134 if (debug)
135 {
136 Pout<< "volPointInterpolation::addSeparated" << endl;
137 }
138
139 auto& pfi = pf.ref();
140 auto& pfbf = pf.boundaryFieldRef();
141
142 const label startOfRequests = UPstream::nRequests();
143
144 forAll(pfbf, patchi)
145 {
146 if (pfbf[patchi].coupled())
147 {
149 (pfbf[patchi]).initSwapAddSeparated
150 (
152 pfi
153 );
154 }
155 }
156
157 // Wait for outstanding requests
158 UPstream::waitRequests(startOfRequests);
159
160 forAll(pfbf, patchi)
161 {
162 if (pfbf[patchi].coupled())
163 {
165 (pfbf[patchi]).swapAddSeparated
166 (
168 pfi
169 );
170 }
171 }
172}
173
174
175template<class Type>
177(
181) const
182{
183 const Field<Type>& pfi = pf.primitiveField();
184
185 // Get face data in flat list
186 const fvMesh& Mesh = mesh();
187 const fvBoundaryMesh& bm = Mesh.boundary();
188
189 auto tboundaryVals = tmp<Field<Type>>::New(Mesh.nBoundaryFaces(), Zero);
190 auto& boundaryVals = tboundaryVals.ref();
191
192 // Do points on 'normal' patches from the surrounding patch faces
193 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
194
195 const primitivePatch& boundary = boundaryPtr_();
196 const labelList& mp = boundary.meshPoints();
197
198 forAll(mp, i)
199 {
200 label pointi = mp[i];
201
202 if (isPatchPoint_[pointi])
203 {
204 const labelList& pFaces = boundary.pointFaces()[i];
205 const scalarList& pWeights = boundaryPointWeights_[i];
206 const Type& val = pfi[pointi];
207 // Face-to-point weights should, in general, have half the weight
208 // of what they actually do in volPointInterpolation since, in
209 // a complete case, a face laying on the opposite side of the
210 // symmetry plane would also contribute to a point laying on
211 // the symmetry plane.
212 // For face-to-point interpolation this is not a problem, but for
213 // the adjoint point-to-face interpolation, the correct value of
214 // the weight should be taken into consideration
215 scalar mod(isSymmetryPoint_[pointi] ? 0.5 : 1);
216
217 forAll(pFaces, j)
218 {
219 if (boundaryIsPatchFace_[pFaces[j]])
220 {
221 boundaryVals[pFaces[j]] += mod*pWeights[j]*val;
222 }
223 }
224 }
225 }
226
227 // Transfer values to face-based sensitivity field
228 for (const label patchi : patchIDs)
229 {
230 label bFacei = bm[patchi].patch().start() - Mesh.nInternalFaces();
231 if (!isA<emptyFvPatch>(bm[patchi]) && !vf[patchi].coupled())
232 {
233 vf[patchi] =
235 (
236 boundaryVals,
237 vf[patchi].size(),
238 bFacei
239 );
240 }
241 }
242}
243
244
245template<class Type>
247(
250) const
251{
252 const primitivePatch& boundary = boundaryPtr_();
253
254 Field<Type>& pfi = pf.primitiveFieldRef();
255
256 // Get face data in flat list
257 tmp<Field<Type>> tboundaryVals(flatBoundaryField(vf));
258 const Field<Type>& boundaryVals = tboundaryVals();
259
260
261 // Do points on 'normal' patches from the surrounding patch faces
262 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
263
264 const labelList& mp = boundary.meshPoints();
265
266 forAll(mp, i)
267 {
268 label pointi = mp[i];
269
270 if (isPatchPoint_[pointi])
271 {
272 const labelList& pFaces = boundary.pointFaces()[i];
273 const scalarList& pWeights = boundaryPointWeights_[i];
274
275 Type& val = pfi[pointi];
276
277 val = Zero;
278 forAll(pFaces, j)
279 {
280 if (boundaryIsPatchFace_[pFaces[j]])
281 {
282 scalar mod(1);
283 if (isSymmetryPoint_[pointi])
284 {
285 const label globalFaceI =
286 mesh().nInternalFaces() + pFaces[j];
287 const label facePatchID =
288 mesh().boundaryMesh().whichPatch(globalFaceI);
289 const polyPatch& pp = mesh().boundaryMesh()[facePatchID];
290 if
291 (
294 )
295 {
296 mod = 0;
297 }
298 //else
299 //{
300 // mod = 2;
301 //}
302 }
303 val += mod*pWeights[j]*boundaryVals[pFaces[j]];
304 }
305 }
306 }
307 }
308
309 // Sum collocated contributions
311
312 // And add separated contributions
313 addSeparated(pf);
314
315 // Optionally normalise
316 /*
317 if (normalisationPtr_)
318 {
319 const scalarField& normalisation = normalisationPtr_();
320 forAll(mp, i)
321 {
322 pfi[mp[i]] *= normalisation[i];
323 }
324 }
325 */
326
327
328 // Push master data to slaves. It is possible (not sure how often) for
329 // a coupled point to have its master on a different patch so
330 // to make sure just push master data to slaves.
331 pushUntransformedData(pfi);
332
333 // Apply displacement constraints
334 const pointConstraints& pcs = pointConstraints::New(pf.mesh());
335
336 pcs.constrain(pf, false);
337}
338
339
340// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
const Mesh & mesh() const noexcept
Return const reference to mesh.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of boundary fields.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
static FOAM_NO_DANGLING_REFERENCE const pointConstraints & New(const pointMesh &mesh, Args &&... args)
const labelList & meshPoints() const
Return labelList of mesh points in patch.
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests).
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
static void waitRequests()
Wait for all requests to finish.
Definition UPstream.H:2497
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
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
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
const mapDistribute & globalCoPointSlavesMap() const
const labelListList & globalCoPointSlaves() const
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
label constructSize() const noexcept
Constructed data size.
Class containing processor-to-processor mapping information.
void reverseDistribute(const label constructSize, List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Reverse distribute data using default commsType.
Application of (multi-)patch point constraints.
void constrain(GeometricField< Type, pointPatchField, pointMesh > &pf, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints) and.
static void syncUntransformedData(const polyMesh &mesh, List< Type > &pointData, const CombineOp &cop)
Helper: sync data on collocated points only.
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const globalMeshData & globalData() const
Return parallel info (demand-driven).
Definition polyMesh.C:1296
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).
label nInternalFaces() const noexcept
Number of internal faces.
A class for managing temporary objects.
Definition tmp.H:75
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 interpolateBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf, GeometricField< Type, pointPatchField, pointMesh > &pf) const
scalarListList boundaryPointWeights_
Per boundary point the weights per pointFaces.
void addSeparated(GeometricField< Type, pointPatchField, pointMesh > &) const
Add separated contributions.
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.
tmp< Field< Type > > flatBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Get boundary field in same order as boundary faces. Field is.
autoPtr< primitivePatch > boundaryPtr_
Boundary addressing.
bitSet boundaryIsPatchFace_
Per boundary face whether is on non-coupled, non-empty patch.
bitSet isSymmetryPoint_
Per mesh(!) point whether is on symmetry plane.
faceListList boundary
bool coupled
dynamicFvMesh & mesh
const dimensionedScalar mp
Proton mass.
Namespace for handling debugging switches.
Definition debug.C:45
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
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.
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
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
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]
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299