Go to the documentation of this file.
40 const scalar coupledPolyPatch::defaultMatchTol_ = 1
e-4;
50 { transformType::UNKNOWN,
"unknown" },
51 { transformType::ROTATIONAL,
"rotational" },
52 { transformType::TRANSLATIONAL,
"translational" },
53 { transformType::COINCIDENTFULLMATCH,
"coincidentFullMatch" },
54 { transformType::NOORDERING,
"noOrdering" },
62 os <<
"v " << pt.
x() <<
' ' << pt.
y() <<
' ' << pt.
z() <<
endl;
94 os <<
"l " << vertI-1 <<
' ' << vertI <<
nl;
111 for (
const face&
f : faces)
115 if (foamToObj.insert(
f[fp], vertI))
125 os <<
' ' << foamToObj[
f[fp]]+1;
127 os <<
' ' << foamToObj[
f[0]]+1 <<
nl;
134 const UList<face>& faces,
146 anchors[facei] =
points[faces[facei][0]];
154 const face&
f = faces[facei];
164 for (label fp2 = 0; fp2 <
f.size(); ++fp2)
166 if (
f[fp1] ==
f[fp2])
191 anchors[facei] =
points[faces[facei][0]];
202 const UList<face>& faces,
212 const point& cc = faceCentres[facei];
214 const face&
f = faces[facei];
218 scalar maxLenSqr = -GREAT;
221 scalar maxCmpt = -GREAT;
226 maxLenSqr =
max(maxLenSqr,
magSqr(pt - cc));
249 scalar minDistSqr = GREAT;
255 if (distSqr < minDistSqr)
257 minDistSqr = distSqr;
262 if (anchorFp == -1 ||
Foam::sqrt(minDistSqr) > tol)
273 if (distSqr == minDistSqr && fp != anchorFp)
276 <<
"Cannot determine unique anchor point on face "
277 << UIndirectList<point>(
points,
f)
279 <<
"Both at index " << anchorFp <<
" and " << fp
280 <<
" the vertices have the same distance "
282 <<
" to the anchor " << anchor
283 <<
". Continuing but results might be wrong."
289 return (
f.size() - anchorFp) %
f.size();
307 Pout<<
"coupledPolyPatch::calcTransformTensors : " <<
name() <<
endl
308 <<
" transform:" << transformTypeNames[
transform] <<
nl
309 <<
" (half)size:" << Cf.size() <<
nl
310 <<
" absTol:" << absTol <<
nl
311 <<
" smallDist min:" <<
min(smallDist) <<
nl
312 <<
" smallDist max:" <<
max(smallDist) <<
nl
313 <<
" sum(mag(nf & nr)):" <<
sum(
mag(nf & nr)) <<
endl;
327 separation_.setSize(0);
330 collocated_.setSize(0);
334 scalar error = absTol*
Foam::sqrt(1.0*Cf.size());
347 && (
sum(
mag(nf & nr)) < Cf.size() - error)
355 separation_.setSize(0);
357 forwardT_.setSize(Cf.size());
358 reverseT_.setSize(Cf.size());
359 collocated_.setSize(Cf.size());
370 Pout<<
" sum(mag(forwardT_ - forwardT_[0])):"
371 <<
sum(
mag(forwardT_ - forwardT_[0]))
375 if (
sum(
mag(forwardT_ - forwardT_[0])) < error)
377 forwardT_.setSize(1);
378 reverseT_.setSize(1);
379 collocated_.setSize(1);
383 Pout<<
" difference in rotation less than"
384 <<
" local tolerance "
385 << error <<
". Assuming uniform rotation." <<
endl;
393 forwardT_.setSize(0);
394 reverseT_.setSize(0);
396 separation_ = Cr - Cf;
398 collocated_.setSize(separation_.size());
406 bool sameSeparation =
true;
407 bool doneWarning =
false;
409 forAll(separation_, facei)
411 scalar smallSqr =
sqr(smallDist[facei]);
413 collocated_[facei] = (
magSqr(separation_[facei]) < smallSqr);
416 if (
magSqr(separation_[facei] - separation_[0]) > smallSqr)
418 sameSeparation =
false;
420 if (!doneWarning &&
debug)
424 Pout<<
" separation " << separation_[facei]
426 <<
" differs from separation[0] " << separation_[0]
427 <<
" by more than local tolerance "
429 <<
". Assuming non-uniform separation." <<
endl;
441 Pout<<
" separation " <<
mag(separation_[0])
442 <<
" less than local tolerance " << smallDist[0]
443 <<
". Assuming zero separation." <<
endl;
446 separation_.setSize(0);
453 Pout<<
" separation " <<
mag(separation_[0])
454 <<
" more than local tolerance " << smallDist[0]
455 <<
". Assuming uniform separation." <<
endl;
458 separation_.setSize(1);
467 Pout<<
" separation_:" << separation_.size() <<
nl
468 <<
" forwardT size:" << forwardT_.size() <<
endl;
481 const polyBoundaryMesh& bm,
482 const word& patchType,
487 matchTolerance_(defaultMatchTol_),
498 const word& patchType
502 matchTolerance_(
dict.getOrDefault(
"matchTolerance", defaultMatchTol_)),
505 transformTypeNames.getOrDefault
522 matchTolerance_(pp.matchTolerance_),
523 transform_(pp.transform_)
534 matchTolerance_(pp.matchTolerance_),
535 transform_(pp.transform_)
549 matchTolerance_(pp.matchTolerance_),
550 transform_(pp.transform_)
563 polyPatch(pp, bm, index, mapAddressing, newStart),
564 matchTolerance_(pp.matchTolerance_),
565 transform_(pp.transform_)
583 os.
writeEntry(
"transform", transformTypeNames[transform_]);
static pointField getAnchorPoints(const UList< face > &, const pointField &, const transformType)
vectorField pointField
pointField is a vectorField.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
A class for handling words, derived from Foam::string.
A class for handling file names.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
static const Enum< transformType > transformTypeNames
List< bool > boolList
A List of bools.
Ostream & endl(Ostream &os)
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
dimensionSet transform(const dimensionSet &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
virtual void write(Ostream &os) const
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Generic templated field type.
A patch is a list of labels that address the faces in the global face list.
static label getRotation(const pointField &points, const face &f, const point &anchor, const scalar tol)
label max(const labelHashSet &set, label maxValue=labelMin)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual void write(Ostream &os) const
OBJstream os(runTime.globalPath()/outputName)
static scalarField calcFaceTol(const UList< face > &faces, const pointField &points, const pointField &faceCentres)
Output to file stream, using an OSstream.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual ~coupledPolyPatch()
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionedScalar e
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
void calcTransformTensors(const vectorField &Cf, const vectorField &Cr, const vectorField &nf, const vectorField &nr, const scalarField &smallDist, const scalar absTol, const transformType=UNKNOWN) const
Ostream & writeEntry(const keyType &key, const T &value)
word name(const expressions::valueTypeCode typeCode)
A List with indirect addressing.
A face is a list of labels corresponding to mesh vertices.
Various functions to operate on Lists.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
const volScalarField & p0
tensor rotationTensor(const vector &n1, const vector &n2)
static void writeOBJ(Ostream &os, const point &pt)
vector point
Point is a vector.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
coupledPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform)
defineTypeNameAndDebug(combustionModel, 0)
#define WarningInFunction
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Smooth ATC in cells next to a set of patches supplied by type.
labelList pointLabels(nPoints, -1)