37 const scalar coupledPolyPatch::defaultMatchTol_ = 1
e-4;
45 "coincidentFullMatch",
58 os <<
"v " << pt.
x() <<
' ' << pt.
y() <<
' ' << pt.
z() <<
endl;
71 writeOBJ(os, points[pointLabels[i]]);
90 os <<
"l " << vertI-1 <<
' ' << vertI <<
nl;
109 const face& f = faces[i];
113 if (foamToObj.insert(f[fp], vertI))
123 os <<
' ' << foamToObj[f[fp]]+1;
125 os <<
' ' << foamToObj[f[0]]+1 <<
nl;
139 if (transform != COINCIDENTFULLMATCH)
144 anchors[faceI] = points[faces[faceI][0]];
152 const face& f = faces[faceI];
158 const point& p1 = points[f[fp1]];
162 for (
label fp2 = 0; fp2 < f.size(); ++fp2)
164 if (f[fp1] == f[fp2])
169 const point& p2 = points[f[fp2]];
189 anchors[faceI] = points[faces[faceI][0]];
210 const point& cc = faceCentres[faceI];
212 const face& f = faces[faceI];
216 scalar maxLenSqr = -GREAT;
219 scalar maxCmpt = -GREAT;
223 const point& pt = points[f[fp]];
224 maxLenSqr =
max(maxLenSqr,
magSqr(pt - cc));
247 scalar minDistSqr = GREAT;
251 scalar distSqr =
magSqr(anchor - points[f[fp]]);
253 if (distSqr < minDistSqr)
255 minDistSqr = distSqr;
260 if (anchorFp == -1 ||
Foam::sqrt(minDistSqr) > tol)
269 scalar distSqr =
magSqr(anchor - points[f[fp]]);
271 if (distSqr == minDistSqr && fp != anchorFp)
275 "label coupledPolyPatch::getRotation\n" 277 " const pointField&,\n" 282 ) <<
"Cannot determine unique anchor point on face " 285 <<
"Both at index " << anchorFp <<
" and " << fp
286 <<
" the vertices have the same distance " 288 <<
" to the anchor " << anchor
289 <<
". Continuing but results might be wrong." 295 return (f.
size() - anchorFp) % f.
size();
313 Pout<<
"coupledPolyPatch::calcTransformTensors : " <<
name() <<
endl 314 <<
" transform:" << transformTypeNames[
transform] <<
nl 315 <<
" (half)size:" << Cf.
size() <<
nl 316 <<
" absTol:" << absTol <<
nl 317 <<
" smallDist min:" <<
min(smallDist) <<
nl 318 <<
" smallDist max:" <<
max(smallDist) <<
nl 319 <<
" sum(mag(nf & nr)):" <<
sum(
mag(nf & nr)) <<
endl;
333 separation_.setSize(0);
336 collocated_.setSize(0);
349 transform == ROTATIONAL
351 transform != TRANSLATIONAL
352 && transform != COINCIDENTFULLMATCH
361 separation_.setSize(0);
363 forwardT_.setSize(Cf.
size());
364 reverseT_.setSize(Cf.
size());
365 collocated_.setSize(Cf.
size());
376 Pout<<
" sum(mag(forwardT_ - forwardT_[0])):" 377 <<
sum(
mag(forwardT_ - forwardT_[0]))
381 if (
sum(
mag(forwardT_ - forwardT_[0])) < error)
383 forwardT_.setSize(1);
384 reverseT_.setSize(1);
385 collocated_.setSize(1);
389 Pout<<
" difference in rotation less than" 390 <<
" local tolerance " 391 << error <<
". Assuming uniform rotation." <<
endl;
399 forwardT_.setSize(0);
400 reverseT_.setSize(0);
402 separation_ = Cr - Cf;
404 collocated_.
setSize(separation_.size());
412 bool sameSeparation =
true;
413 bool doneWarning =
false;
415 forAll(separation_, facei)
417 scalar smallSqr =
sqr(smallDist[facei]);
419 collocated_[facei] = (
magSqr(separation_[facei]) < smallSqr);
422 if (
magSqr(separation_[facei] - separation_[0]) > smallSqr)
424 sameSeparation =
false;
426 if (!doneWarning && debug)
430 Pout<<
" separation " << separation_[facei]
432 <<
" differs from separation[0] " << separation_[0]
433 <<
" by more than local tolerance " 435 <<
". Assuming non-uniform separation." <<
endl;
447 Pout<<
" separation " <<
mag(separation_[0])
448 <<
" less than local tolerance " << smallDist[0]
449 <<
". Assuming zero separation." <<
endl;
452 separation_.setSize(0);
459 Pout<<
" separation " <<
mag(separation_[0])
460 <<
" more than local tolerance " << smallDist[0]
461 <<
". Assuming uniform separation." <<
endl;
464 separation_.setSize(1);
473 Pout<<
" separation_:" << separation_.size() <<
nl 474 <<
" forwardT size:" << forwardT_.size() <<
endl;
488 const word& patchType,
492 polyPatch(name, size, start, index, bm, patchType),
493 matchTolerance_(defaultMatchTol_),
494 transform_(transform)
504 const word& patchType
507 polyPatch(name, dict, index, bm, patchType),
508 matchTolerance_(dict.
lookupOrDefault(
"matchTolerance", defaultMatchTol_)),
511 dict.
found(
"transform")
512 ? transformTypeNames.read(dict.
lookup(
"transform"))
525 matchTolerance_(pp.matchTolerance_),
526 transform_(pp.transform_)
539 polyPatch(pp, bm, index, newSize, newStart),
540 matchTolerance_(pp.matchTolerance_),
541 transform_(pp.transform_)
554 polyPatch(pp, bm, index, mapAddressing, newStart),
555 matchTolerance_(pp.matchTolerance_),
556 transform_(pp.transform_)
575 os.
writeKeyword(
"transform") << transformTypeNames[transform_]
dimensionedScalar sqrt(const dimensionedScalar &ds)
static void writeOBJ(Ostream &os, const point &pt)
Write point in OBJ format.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
const pointField & points
dimensioned< scalar > mag(const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
virtual ~coupledPolyPatch()
Destructor.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
A class for handling words, derived from string.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
coupledPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform)
Construct from components.
void size(const label)
Override size to be inconsistent with allocated storage.
Various functions to operate on Lists.
A list of keyword definitions, which are a keyword followed by any number of values (e...
A patch is a list of labels that address the faces in the global face list.
A face is a list of labels corresponding to mesh vertices.
const double e
Elementary charge.
void setSize(const label)
Reset size of List.
Ostream & endl(Ostream &os)
Add newline and flush stream.
static pointField getAnchorPoints(const UList< face > &, const pointField &, const transformType)
Get a unique anchor point for all faces.
#define WarningIn(functionName)
Report a warning using Foam::Warning.
A List with indirect addressing.
dimensionSet transform(const dimensionSet &)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
label size() const
Return the number of elements in the UList.
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
Calculate the transformation tensors.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
static label getRotation(const pointField &points, const face &f, const point &anchor, const scalar tol)
Get the number of vertices face f needs to be rotated such that.
Initialise the NamedEnum HashTable from the static list of names.
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Class to handle errors and exceptions in a simple, consistent stream-based manner.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
A class for handling file names.
static scalarField calcFaceTol(const UList< face > &faces, const pointField &points, const pointField &faceCentres)
Calculate typical tolerance per face. Is currently max distance.
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
3D tensor transformation operations.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
static const NamedEnum< transformType, 5 > transformTypeNames
tensor rotationTensor(const vector &n1, const vector &n2)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout,"Pout")
List< bool > boolList
Bool container classes.