Go to the documentation of this file.
50 Pout<<
"volPointInterpolation::calcBoundaryAddressing() : "
51 <<
"constructing boundary addressing"
62 mesh().nFaces()-
mesh().nInternalFaces(),
63 mesh().nInternalFaces()
88 !isA<emptyPolyPatch>(pp)
125 Pout<<
"volPointInterpolation::calcBoundaryAddressing():"
126 <<
" added dangling mesh point:" << pointI
132 label nPatchFace = 0;
140 label nPatchPoint = 0;
150 <<
" of which on proper patch:" << nPatchFace <<
nl
152 <<
" of which on proper patch:" << nPatchPoint <<
endl;
161 Pout<<
"volPointInterpolation::makeInternalWeights() : "
162 <<
"constructing weighting factors for internal and non-coupled"
163 <<
" points." <<
endl;
180 const labelList& pcp = pointCells[pointi];
188 1.0/
mag(
points[pointi] - cellCentres[pcp[pointCelli]]);
190 sumWeights[pointi] += pw[pointCelli];
201 Pout<<
"volPointInterpolation::makeBoundaryWeights() : "
202 <<
"constructing weighting factors for boundary points." <<
endl;
224 sumWeights[pointI] = 0.0;
232 pw[i] = 1.0/
mag(
points[pointI] - faceCentres[faceI]);
233 sumWeights[pointI] += pw[i];
249 Pout<<
"volPointInterpolation::makeWeights() : "
250 <<
"constructing weighting factors"
263 "volPointSumWeights",
320 pw[i] /= sumWeights[pointI];
335 pw[i] /= sumWeights[pointI];
342 Pout<<
"volPointInterpolation::makeWeights() : "
343 <<
"finished constructing weighting factors"
387 interpolateInternalField(vf, pf);
390 interpolateBoundaryField(vf, pf);
virtual const pointField & points() const
Return raw points.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
#define forAll(list, i)
Loop across all elements in list.
boolList boundaryIsPatchFace_
Per boundary face whether is on non-coupled patch.
void makeBoundaryWeights(scalarField &sumWeights)
Make weights for points on uncoupled patches.
Template functions to aid in the implementation of demand driven data.
A List obtained as a section of another List.
volPointInterpolation(const volPointInterpolation &)
Disallow default bitwise copy construct.
autoPtr< primitivePatch > boundaryPtr_
Boundary addressing.
faceListList boundary(nPatches)
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Ostream & endl(Ostream &os)
Add newline and flush stream.
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
dimensioned< scalar > mag(const dimensioned< Type > &)
const fileName & instance() const
void updateMesh(const mapPolyMesh &)
Update mesh topology using the morph engine.
static const pointMesh & New(const polyMesh &mesh)
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("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]=cellShape(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]
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
void makeWeights()
Construct all point weighting factors.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
const fvMesh & mesh() const
A patch is a list of labels that address the faces in the global face list.
void addSeparated(GeometricField< Type, pointPatchField, pointMesh > &) const
Add separated contributions.
void calcBoundaryAddressing()
Construct addressing over all boundary faces.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label nInternalFaces() const
Mesh data needed to do the Finite Volume discretisation.
~volPointInterpolation()
Destructor.
label start() const
Return start label of this patch in the polyMesh face list.
Application of (multi-)patch point contraints.
void setSize(const label)
Reset size of List.
const vectorField & cellCentres() const
prefixOSstream Pout(cout, "Pout")
static void syncUntransformedData(const polyMesh &mesh, List< Type > &pointData, const CombineOp &cop)
Helper: sync data on collocated points only.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
bool movePoints()
Correct weighting factors for moving mesh.
const vectorField & faceCentres() const
void clear()
Clear the list, i.e. set size to zero.
const labelListList & pointCells() const
void interpolateDisplacement(const volVectorField &, pointVectorField &) const
Interpolate from volField to pointField.
void makeInternalWeights(scalarField &sumWeights)
Make weights for internal and coupled-only boundarypoints.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
A face is a list of labels corresponding to mesh vertices.
void size(const label)
Override size to be inconsistent with allocated storage.
scalarListList boundaryPointWeights_
Per boundary point the weights per pointFaces.
scalarListList pointWeights_
Interpolation scheme weighting factor array.
boolList isPatchPoint_
Per mesh(!) point whether is on non-coupled patch (on any.
Generic GeometricField class.
defineTypeNameAndDebug(combustionModel, 0)
Interpolate from cell centres to points (vertices) using inverse distance weighting.
void pushUntransformedData(List< Type > &) const
Helper: push master point data to collocated points.
A list of faces which address into the list of points.