Go to the documentation of this file.
55 const vectorField n((faceCentres - rotationCentre_) ^ rotationAxis_);
63 Info<<
"findFaceMaxRadius(const pointField&) : patch: " <<
name() <<
nl
64 <<
" rotFace = " << faceI <<
nl
65 <<
" point = " << faceCentres[faceI] <<
nl
66 <<
" distance = " <<
Foam::sqrt(magRadSqr[faceI])
82 half0Areas[facei] = half0[facei].normal(half0.
points());
89 half1Areas[facei] = half1[facei].normal(half1.
points());
103 Pout<<
"calcTransforms() : patch: " <<
name() <<
nl
124 <<
"Patch " <<
name()
125 <<
" has transform type " << transformTypeNames[
transform()]
126 <<
", neighbour patch " << neighbPatchName()
127 <<
" has transform type "
128 << neighbPatch().transformTypeNames[neighbPatch().transform()]
141 if (rotationAngleDefined_)
143 const tensor T(rotationAxis_*rotationAxis_);
147 0, -rotationAxis_.z(), rotationAxis_.y(),
148 rotationAxis_.z(), 0, -rotationAxis_.x(),
149 -rotationAxis_.y(), rotationAxis_.x(), 0
156 +
sin(rotationAngle_)*S
163 +
sin(-rotationAngle_)*S
168 const vector transformedAreaPos =
gSum(half1Areas & revTPos);
169 const vector transformedAreaNeg =
gSum(half1Areas & revTNeg);
171 const scalar magArea0 =
mag(area0) + ROOTVSMALL;
175 const scalar errorPos =
mag(transformedAreaPos + area0);
176 const scalar errorNeg =
mag(transformedAreaNeg + area0);
178 const scalar normErrorPos = errorPos/magArea0;
179 const scalar normErrorNeg = errorNeg/magArea0;
181 if (errorPos > errorNeg && normErrorNeg < matchTolerance())
184 rotationAngle_ *= -1;
191 const scalar areaError =
min(normErrorPos, normErrorNeg);
193 if (areaError > matchTolerance())
196 <<
"Patch areas are not consistent within "
197 << 100*matchTolerance()
198 <<
" % indicating a possible error in the specified "
199 <<
"angle of rotation" <<
nl
200 <<
" owner patch : " <<
name() <<
nl
201 <<
" neighbour patch : " << neighbPatch().name()
205 <<
" area error : " << 100*areaError <<
" %"
206 <<
" match tolerance : " << matchTolerance()
212 scalar theta =
radToDeg(rotationAngle_);
214 Pout<<
"cyclicAMIPolyPatch::calcTransforms: patch:"
216 <<
" Specified rotation:"
217 <<
" swept angle: " << theta <<
" [deg]"
218 <<
" reverse transform: " << revT
226 if (half0Ctrs.size())
228 n0 = findFaceNormalMaxRadius(half0Ctrs);
230 if (half1Ctrs.size())
232 n1 = -findFaceNormalMaxRadius(half1Ctrs);
238 n0 /=
mag(n0) + VSMALL;
239 n1 /=
mag(n1) + VSMALL;
246 (n0 ^ rotationAxis_),
252 (-n1 ^ rotationAxis_),
261 Pout<<
"cyclicAMIPolyPatch::calcTransforms: patch:"
263 <<
" Specified rotation:"
264 <<
" n0:" << n0 <<
" n1:" << n1
265 <<
" swept angle: " << theta <<
" [deg]"
266 <<
" reverse transform: " << revT
282 Pout<<
"cyclicAMIPolyPatch::calcTransforms : patch:" <<
name()
283 <<
" Specified translation : " << separationVector_
303 <<
" Assuming cyclic AMI pairs are colocated" <<
endl;
318 <<
" forwardT = " << forwardT() <<
nl
319 <<
" reverseT = " << reverseT() <<
nl
320 <<
" separation = " << separation() <<
nl
321 <<
" collocated = " << collocated() <<
nl <<
endl;
341 neighbPatch().meshPoints()
352 transformPosition(nbrPoints);
384 AMILowWeightCorrection_,
391 Pout<<
"cyclicAMIPolyPatch : " <<
name()
392 <<
" constructed AMI with " <<
nl
393 <<
" " <<
"srcAddress:" << AMIPtr_().srcAddress().size()
395 <<
" " <<
"tgAddress :" << AMIPtr_().tgtAddress().size()
419 neighbPatch().faceCentres(),
420 neighbPatch().faceAreas(),
421 neighbPatch().faceCellCentres()
485 const word& patchType,
494 rotationAngleDefined_(
false),
499 AMIRequireMatch_(
true),
500 AMILowWeightCorrection_(-1.0),
515 const word& patchType
524 rotationAngleDefined_(
false),
529 AMIRequireMatch_(
true),
539 ) <<
"No \"neighbourPatch\" or \"coupleGroup\" provided."
543 if (nbrPatchName_ ==
name)
548 ) <<
"Neighbour patch name " << nbrPatchName_
549 <<
" cannot be the same as this patch " <<
name
558 dict.
lookup(
"rotationCentre") >> rotationCentre_;
561 rotationAngleDefined_ =
true;
562 rotationAngle_ =
degToRad(rotationAngle_);
566 Info<<
"rotationAngle: " << rotationAngle_ <<
" [rad]"
571 scalar magRot =
mag(rotationAxis_);
577 ) <<
"Illegal rotationAxis " << rotationAxis_ <<
endl
578 <<
"Please supply a non-zero vector."
581 rotationAxis_ /= magRot;
587 dict.
lookup(
"separationVector") >> separationVector_;
634 const label newStart,
635 const word& nbrPatchName
639 nbrPatchName_(nbrPatchName),
654 if (nbrPatchName_ ==
name())
657 <<
"Neighbour patch name " << nbrPatchName_
658 <<
" cannot be the same as this patch " <<
name()
704 if (nbrPatchID_ == -1)
708 if (nbrPatchID_ == -1)
711 <<
"Illegal neighbourPatch name " << neighbPatchName()
712 <<
nl <<
"Valid patch names are "
719 refCast<const cyclicAMIPolyPatch>
727 <<
"Patch " <<
name()
728 <<
" specifies neighbour patch " << neighbPatchName()
729 <<
nl <<
" but that in return specifies "
740 return index() < neighbPatchID();
747 return refCast<const cyclicAMIPolyPatch>(pp);
754 const word surfType(surfDict_.lookupOrDefault<
word>(
"type",
"none"));
756 if (!surfPtr_.valid() && owner() && surfType !=
"none")
758 word surfName(surfDict_.lookupOrDefault(
"name",
name()));
788 <<
"AMI interpolator only available to owner patch"
792 if (!AMIPtr_.valid())
805 return AMI().applyLowWeightCorrection();
809 return neighbPatch().AMI().applyLowWeightCorrection();
828 else if (separated())
859 forwardT().size() == 1
873 else if (separated())
877 separation().size() == 1
879 : separation()[faceI]
897 reverseT().size() == 1
911 else if (separated())
915 separation().size() == 1
917 : separation()[faceI]
935 reverseT().size() == 1
1001 reverseTransformPosition(prt, faceI);
1004 reverseTransformDirection(nrt, faceI);
1006 label nbrFaceI = -1;
1010 nbrFaceI = AMI().tgtPointFace
1021 nbrFaceI = neighbPatch().AMI().srcPointFace
1043 if (!nbrPatchName_.empty())
1048 coupleGroup_.write(os);
1059 if (rotationAngleDefined_)
1069 os.
writeKeyword(
"separationVector") << separationVector_
1089 if (AMILowWeightCorrection_ > 0)
1091 os.
writeKeyword(
"lowWeightCorrection") << AMILowWeightCorrection_
1095 if (!surfDict_.empty())
interpolationMethod
Enumeration specifying interpolation method.
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
const bMesh & mesh() const
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
points setSize(newPointi)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
static bool valid(char)
Is this character valid for a word.
label pointFace(const label faceI, const vector &n, point &p) const
Return face index on neighbour patch which shares point p.
virtual const vectorField & separation() const
If the planes are separated the separation vector.
label findMax(const ListType &, const label start=0)
Find index of max element (and larger than given element).
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
const Field< PointType > & points() const
Return reference to global points.
virtual bool owner() const
Does this side own the patch?
A class for handling words, derived from string.
A class for handling file names.
virtual void movePoints(PstreamBuffers &, const pointField &p)
Correct patches after moving points.
#define forAll(list, i)
Loop across all elements in list.
word nbrPatchName_
Name of other half.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
cyclicAMIPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform=UNKNOWN)
Construct from (base couped patch) components.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
point rotationCentre_
Point on axis of rotation for rotational cyclics.
dimensionedScalar sin(const dimensionedScalar &ds)
A List obtained as a section of another List.
virtual const tensorField & forwardT() const
Return face transformation tensor.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
virtual void clearGeom()
Clear geometry.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
Combination-Reduction operation for a parallel run.
const dictionary surfDict_
Dictionary used during projection surface construction.
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
const autoPtr< searchableSurface > & surfPtr() const
Return a reference to the projection surface.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const coupleGroupIdentifier coupleGroup_
Optional patchGroup to find neighbPatch.
Ostream & endl(Ostream &os)
Add newline and flush stream.
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
Type gSum(const FieldField< Field, Type > &f)
virtual void reverseTransformPosition(point &l, const label faceI) const
Transform a patch-based position from this side to nbr side.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Base class for Arbitrary Mesh Interface (AMI) methods.
Mesh consisting of general polyhedral cells.
dimensionSet transform(const dimensionSet &)
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
vector findFaceNormalMaxRadius(const pointField &faceCentres) const
Return normal of face at max distance from rotation axis.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
virtual void transformPosition(pointField &) const
Transform patch-based positions from nbr side to this side.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
vector separationVector_
Translation vector.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Pre-declare SubField and related Field type.
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
A patch is a list of labels that address the faces in the global face list.
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual void movePoints(const Field< PointType > &)
Correct patch after moving points.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
virtual void calcTransforms()
Recalculate the transformation tensors.
vector rotationAxis_
Axis of rotation for rotational cyclics.
scalar rotationAngle_
Rotation angle.
virtual label neighbPatchID() const
Neighbour patch ID.
A list of keyword definitions, which are a keyword followed by any number of values (e....
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
Macros for easy insertion into run-time selection tables.
bool rotationAngleDefined_
Flag to show whether the rotation angle is defined.
errorManip< error > abort(error &err)
List< bool > boolList
Bool container classes.
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
virtual void reverseTransformDirection(vector &d, const label faceI) const
Transform a patch-based direction from this side to nbr side.
errorManipArg< error, int > exit(error &err, const int errNo=1)
virtual void resetAMI(const AMIPatchToPatchInterpolation::interpolationMethod &AMIMethod=AMIPatchToPatchInterpolation::imFaceAreaWeight) const
Reset the AMI interpolator.
void setSize(const label)
Reset size of List.
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
fileName path() const
Return path.
prefixOSstream Pout(cout, "Pout")
const vectorField::subField faceCentres() const
Return face centres.
bool applyLowWeightCorrection() const
Return true if applying the low weight correction.
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
virtual void clearGeom()
Clear geometry.
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)
dimensionedScalar acos(const dimensionedScalar &ds)
const word & constant() const
Return constant name.
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
const word & neighbPatchName() const
Neighbour patch name.
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
static const word null
An empty word.
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
const Time & time() const
Return the top-level database.
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
label findPatchID(const polyPatchList &, const word &) const
Get index of polypatch by name.
virtual const boolList & collocated() const
Are faces collocated. Either size 0,1 or length of patch.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
const scalar AMILowWeightCorrection_
Low weight correction threshold for AMI.
bool AMIRequireMatch_
Flag to indicate that patches should match/overlap.
defineTypeNameAndDebug(combustionModel, 0)
const bool AMIReverse_
Flag to indicate that slave patch should be reversed for AMI.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const word & name() const
Return name.
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
virtual ~cyclicAMIPolyPatch()
Destructor.
word name(const complex &)
Return a string representation of a complex.
A list of faces which address into the list of points.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
dimensionedScalar cos(const dimensionedScalar &ds)
Cyclic patch for Arbitrary Mesh Interface (AMI)
Tensor< Cmpt > T() const
Transpose.