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 -------------------------------------------------------------------------------
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 "volFields.H"
31 #include "pointFields.H"
32 #include "emptyFvPatch.H"
33 #include "coupledPointPatchField.H"
34 #include "pointConstraints.H"
36 #include "symmetryPlanePolyPatch.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 template<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 
83 template<class Type>
85 (
86  const GeometricField<Type, fvPatchField, volMesh>& vf
87 ) const
88 {
89  const fvMesh& mesh = vf.mesh();
90  const fvBoundaryMesh& bm = mesh.boundary();
91 
92  tmp<Field<Type>> tboundaryVals
93  (
94  new Field<Type>(mesh.nBoundaryFaces())
95  );
96  Field<Type>& boundaryVals = tboundaryVals.ref();
97 
98  forAll(vf.boundaryField(), patchi)
99  {
100  label bFacei = bm[patchi].patch().start() - mesh.nInternalFaces();
101 
102  if
103  (
104  !isA<emptyFvPatch>(bm[patchi])
105  && !vf.boundaryField()[patchi].coupled()
106  )
107  {
108  SubList<Type>
109  (
110  boundaryVals,
111  vf.boundaryField()[patchi].size(),
112  bFacei
113  ) = vf.boundaryField()[patchi];
114  }
115  else
116  {
117  const polyPatch& pp = bm[patchi].patch();
118 
119  forAll(pp, i)
120  {
121  boundaryVals[bFacei++] = Zero;
122  }
123  }
124  }
126  return tboundaryVals;
127 }
128 
129 
130 template<class Type>
132 (
134 ) const
135 {
136  if (debug)
137  {
138  Pout<< "volPointInterpolation::addSeparated" << endl;
139  }
140 
141  auto& pfi = pf.ref();
142  auto& pfbf = pf.boundaryFieldRef();
143 
144  const label startOfRequests = UPstream::nRequests();
145 
146  forAll(pfbf, patchi)
147  {
148  if (pfbf[patchi].coupled())
149  {
150  refCast<coupledPointPatchField<Type>>
151  (pfbf[patchi]).initSwapAddSeparated
152  (
154  pfi
155  );
156  }
157  }
158 
159  // Wait for outstanding requests
160  UPstream::waitRequests(startOfRequests);
161 
162  forAll(pfbf, patchi)
163  {
164  if (pfbf[patchi].coupled())
165  {
166  refCast<coupledPointPatchField<Type>>
167  (pfbf[patchi]).swapAddSeparated
168  (
170  pfi
171  );
172  }
173  }
174 }
175 
176 
177 template<class Type>
179 (
180  const GeometricField<Type, pointPatchField, pointMesh>& pf,
182  const labelHashSet& patchIDs
183 ) const
184 {
185  const Field<Type>& pfi = pf.primitiveField();
186 
187  // Get face data in flat list
188  const fvMesh& Mesh = mesh();
189  const fvBoundaryMesh& bm = Mesh.boundary();
190 
191  tmp<Field<Type>> tboundaryVals
192  (
193  new Field<Type>(Mesh.nBoundaryFaces(), Zero)
194  );
195  Field<Type>& boundaryVals = tboundaryVals.ref();
196 
197  // Do points on 'normal' patches from the surrounding patch faces
198  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
199 
200  const primitivePatch& boundary = boundaryPtr_();
201  const labelList& mp = boundary.meshPoints();
202 
203  forAll(mp, i)
204  {
205  label pointi = mp[i];
206 
207  if (isPatchPoint_[pointi])
208  {
209  const labelList& pFaces = boundary.pointFaces()[i];
210  const scalarList& pWeights = boundaryPointWeights_[i];
211  const Type& val = pfi[pointi];
212  // Face-to-point weights should, in general, have half the weight
213  // of what they actually do in volPointInterpolation since, in
214  // a complete case, a face laying on the opposite side of the
215  // symmetry plane would also contribute to a point laying on
216  // the symmetry plane.
217  // For face-to-point interpolation this is not a problem, but for
218  // the adjoint point-to-face interpolation, the correct value of
219  // the weight should be taken into consideration
220  scalar mod(isSymmetryPoint_[pointi] ? 0.5 : 1);
221 
222  forAll(pFaces, j)
223  {
224  if (boundaryIsPatchFace_[pFaces[j]])
225  {
226  boundaryVals[pFaces[j]] += mod*pWeights[j]*val;
227  }
228  }
229  }
230  }
231 
232  // Transfer values to face-based sensitivity field
233  for (const label patchi : patchIDs)
234  {
235  label bFacei = bm[patchi].patch().start() - Mesh.nInternalFaces();
236  if (!isA<emptyFvPatch>(bm[patchi]) && !vf[patchi].coupled())
237  {
238  vf[patchi] =
239  SubList<Type>
240  (
241  boundaryVals,
242  vf[patchi].size(),
243  bFacei
244  );
245  }
246  }
247 }
248 
249 
250 template<class Type>
252 (
253  const GeometricField<Type, fvPatchField, volMesh>& vf,
254  GeometricField<Type, pointPatchField, pointMesh>& pf
255 ) const
256 {
257  const primitivePatch& boundary = boundaryPtr_();
258 
259  Field<Type>& pfi = pf.primitiveFieldRef();
260 
261  // Get face data in flat list
262  tmp<Field<Type>> tboundaryVals(flatBoundaryField(vf));
263  const Field<Type>& boundaryVals = tboundaryVals();
264 
265 
266  // Do points on 'normal' patches from the surrounding patch faces
267  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
268 
269  const labelList& mp = boundary.meshPoints();
270 
271  forAll(mp, i)
272  {
273  label pointi = mp[i];
274 
275  if (isPatchPoint_[pointi])
276  {
277  const labelList& pFaces = boundary.pointFaces()[i];
278  const scalarList& pWeights = boundaryPointWeights_[i];
279 
280  Type& val = pfi[pointi];
281 
282  val = Zero;
283  forAll(pFaces, j)
284  {
285  if (boundaryIsPatchFace_[pFaces[j]])
286  {
287  scalar mod(1);
288  if (isSymmetryPoint_[pointi])
289  {
290  const label globalFaceI =
291  mesh().nInternalFaces() + pFaces[j];
292  const label facePatchID =
293  mesh().boundaryMesh().whichPatch(globalFaceI);
294  const polyPatch& pp = mesh().boundaryMesh()[facePatchID];
295  if
296  (
297  isA<symmetryPolyPatch>(pp)
298  || isA<symmetryPlanePolyPatch>(pp)
299  )
300  {
301  mod = 0;
302  }
303  //else
304  //{
305  // mod = 2;
306  //}
307  }
308  val += mod*pWeights[j]*boundaryVals[pFaces[j]];
309  }
310  }
311  }
312  }
313 
314  // Sum collocated contributions
315  pointConstraints::syncUntransformedData(mesh(), pfi, plusEqOp<Type>());
316 
317  // And add separated contributions
318  addSeparated(pf);
319 
320  // Optionally normalise
321  /*
322  if (normalisationPtr_)
323  {
324  const scalarField& normalisation = normalisationPtr_();
325  forAll(mp, i)
326  {
327  pfi[mp[i]] *= normalisation[i];
328  }
329  }
330  */
331 
332 
333  // Push master data to slaves. It is possible (not sure how often) for
334  // a coupled point to have its master on a different patch so
335  // to make sure just push master data to slaves.
336  pushUntransformedData(pfi);
337 
338  // Apply displacement constraints
339  const pointConstraints& pcs = pointConstraints::New(pf.mesh());
340 
341  pcs.constrain(pf, false);
342 }
343 
344 
345 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
faceListList boundary
const labelList patchIDs(pbm.indices(polyPatchNames, true))
const mapDistribute & globalCoPointSlavesMap() const
static const pointConstraints & New(const pointMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
label constructSize() const noexcept
Constructed data size.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of boundary fields.
static void waitRequests()
Wait for all requests to finish.
Definition: UPstream.H:1538
const labelList & meshPoints() const
Return labelList of mesh points in patch.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
A list of faces which address into the list of points.
void addSeparated(GeometricField< Type, pointPatchField, pointMesh > &) const
Add separated contributions.
void interpolateBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf, GeometricField< Type, pointPatchField, pointMesh > &pf) const
dynamicFvMesh & mesh
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
static void syncUntransformedData(const polyMesh &mesh, List< Type > &pointData, const CombineOp &cop)
Helper: sync data on collocated points only.
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
label nInternalFaces() const noexcept
Number of internal faces.
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
void reverseDistribute(const label constructSize, List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Reverse distribute data using default commsType.
int debug
Static debugging option.
void pushUntransformedData(List< Type > &) const
Helper: push master point data to collocated points.
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 indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
Class containing processor-to-processor mapping information.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
"nonBlocking" : (MPI_Isend, MPI_Irecv)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
tmp< Field< Type > > flatBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Get boundary field in same order as boundary faces. Field is.
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
bool coupled
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)
void constrain(GeometricField< Type, pointPatchField, pointMesh > &pf, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints) and.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const dimensionedScalar mp
Proton mass.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const labelListList & globalCoPointSlaves() const
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127