36 void Foam::motionSmootherAlgo::checkConstraints
38 GeometricField<Type, pointPatchField, pointMesh>& pf
41 typedef GeometricField<Type, pointPatchField, pointMesh> FldType;
43 const polyMesh& mesh = pf.mesh();
45 const polyBoundaryMesh& bm = mesh.boundaryMesh();
49 label nPatchPatchPoints = 0;
53 if (!isA<emptyPolyPatch>(bm[patchi]))
55 nPatchPatchPoints += bm[
patchi].boundaryPoints().size();
60 typename FldType::Boundary& bFld = pf.boundaryField();
78 Field<Type> boundaryPointValues(nPatchPatchPoints);
79 nPatchPatchPoints = 0;
83 if (!isA<emptyPolyPatch>(bm[patchi]))
90 label ppp = meshPoints[bp[pointi]];
91 boundaryPointValues[nPatchPatchPoints++] = pf[ppp];
104 nPatchPatchPoints = 0;
108 if (!isA<emptyPolyPatch>(bm[patchi]))
115 label ppp = meshPoints[bp[pointi]];
117 const Type& savedVal = boundaryPointValues[nPatchPatchPoints++];
119 if (savedVal != pf[ppp])
122 <<
"Patch fields are not consistent on mesh point " 123 << ppp <<
" coordinate " << mesh.points()[ppp]
124 <<
" at patch " << bm[
patchi].name() <<
'.' 126 <<
"Reverse evaluation gives value " << savedVal
127 <<
" , forward evaluation gives value " << pf[ppp]
138 Foam::motionSmootherAlgo::avg
148 "avg("+fld.
name()+
')',
170 if (isMasterEdge_.
get(edgeI) == 1)
172 const edge&
e = edges[edgeI];
173 const scalar w = edgeWeight[edgeI];
175 res[e[0]] += w*fld[e[1]];
176 sumWeight[e[0]] += w;
178 res[e[1]] += w*fld[e[0]];
179 sumWeight[e[1]] += w;
208 if (
mag(sumWeight[pointi]) < vSmall)
211 res[pointi] = fld[pointi];
215 res[pointi] /= sumWeight[pointi];
239 if (isInternalPoint(pointi))
241 newFld[pointi] = 0.5*fld[pointi] + 0.5*avgFld[pointi];
250 template<
class Type,
class CombineOp>
251 void Foam::motionSmootherAlgo::testSyncField
254 const CombineOp& cop,
261 Pout<<
"testSyncField : testing synchronisation of Field<Type>." 277 if (
mag(syncedFld[i] - fld[i]) > maxMag)
280 <<
"On element " << i <<
" value:" << fld[i]
281 <<
" synchronised value:" << syncedFld[i]
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
unsigned int get(const label) const
Get value at index I.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Ostream & endl(Ostream &os)
Add newline and flush stream.
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Generic dimensioned Type class.
const dimensionSet & dimensions() const
Return dimensions.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
List< label > labelList
A List of labels.
errorManip< error > abort(error &err)
void smooth(const GeometricField< Type, pointPatchField, pointMesh > &fld, const scalarField &edgeWeight, GeometricField< Type, pointPatchField, pointMesh > &newFld) const
Fully explicit smoothing of fields (not positions)
const Mesh & mesh() const
Return mesh.
static pointConstraints & New(pointMesh &mesh)
Internal & ref()
Return a reference to the dimensioned internal field.
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
prefixOSstream Pout(cout, "Pout")
dimensioned< scalar > mag(const dimensioned< Type > &)
const doubleScalar e
Elementary charge.
Mesh consisting of general polyhedral cells.
A class for managing temporary objects.