53 Pout<<
"volPointInterpolationAdjoint::calcBoundaryAddressing() : "
54 <<
"constructing boundary addressing"
78 const polyPatch&
pp =
pbm[patchi];
83 && !magSf.boundaryField()[patchi].coupled()
90 for (
const auto&
f :
pp.faces())
92 UIndirectList<bool>(isPatchPoint,
f) =
true;
97 UIndirectList<bool>(isSymmetryPoint,
pp.meshPoints()) =
true;
126 <<
" of which on proper patch:" << nPatchFace <<
nl
128 <<
" of which on proper patch:" << nPatchPoint <<
endl;
140 Pout<<
"volPointInterpolationAdjoint::makeBoundaryWeights() : "
141 <<
"constructing weighting factors for boundary points." <<
endl;
149 boundaryPointWeights_.clear();
150 boundaryPointWeights_.setSize(
boundary.meshPoints().
size());
154 label pointi =
boundary.meshPoints()[i];
156 if (isPatchPoint_[pointi])
161 pw.setSize(
pFaces.size());
163 sumWeights[pointi] = 0.0;
167 if (boundaryIsPatchFace_[
pFaces[i]])
171 pw[i] = 1.0/
mag(
points[pointi] - faceCentres[facei]);
172 sumWeights[pointi] += pw[i];
188 Pout<<
"volPointInterpolationAdjoint::makeWeights() : "
189 <<
"constructing weighting factors"
196 calcBoundaryAddressing();
206 "volPointSumWeights",
218 makeBoundaryWeights(sumWeights);
238 pushUntransformedData(sumWeights);
243 const label pointi =
mp[i];
249 pw[i] /= sumWeights[pointi];
255 Pout<<
"volPointInterpolationAdjoint::makeWeights() : "
256 <<
"finished constructing weighting factors"
314 label pointi = mp[i];
316 if (isPatchPoint_[pointi])
319 const scalarList& pWeights = boundaryPointWeights_[i];
320 const vector& val = pf[pointi];
329 scalar mod(isSymmetryPoint_[pointi] ? 0.5 : 1);
333 if (boundaryIsPatchFace_[
pFaces[j]])
335 boundaryVals[
pFaces[j]] += mod*pWeights[j]*val;
342 label nPassedFaces(0);
349 nPassedFaces +=
patch.size();
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
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,...
const fileName & instance() const noexcept
Read access to instance path component.
void setSize(label n)
Alias for resize().
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
const fvMesh & mesh() const noexcept
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).
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
void set(const bitSet &bitset)
Set specified bits from another bitset.
Mesh data needed to do the Finite Volume discretisation.
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
A range or interval of labels defined by a start and a size.
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.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
virtual const pointField & points() const
Return raw points.
A patch is a list of labels that address the faces in the global face list.
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.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
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.
~volPointInterpolationAdjoint()
Destructor.
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.
const dimensionedScalar mp
Proton mass.
Namespace for handling debugging switches.
const std::string patch
OpenFOAM patch number as a std::string.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
List< label > labelList
A List of labels.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
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.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
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.
static constexpr const zero Zero
Global zero (0).
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.
List< scalar > scalarList
List of scalar.
constexpr char nl
The newline '\n' character (0x0a).
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.