62void Foam::cyclicFaPatch::calcTransforms()
64 const label sizeby2 = this->size()/2;
69 for (label edgei=0; edgei < sizeby2; ++edgei)
72 half1Ctrs[edgei] = this->
edgeCentres()[edgei + sizeby2];
78 const vectorField eN(edgeNormals()*magEdgeLengths());
80 scalar maxMatchError = 0;
83 for (label edgei = 0; edgei < sizeby2; ++edgei)
85 half0Normals[edgei] = eN[edgei];
86 half1Normals[edgei] = eN[edgei + sizeby2];
88 scalar magLe =
mag(half0Normals[edgei]);
89 scalar nbrMagLe =
mag(half1Normals[edgei]);
90 scalar avLe = 0.5*(magLe + nbrMagLe);
92 if (magLe < ROOTVSMALL && nbrMagLe < ROOTVSMALL)
97 half0Normals[edgei] =
point(1, 0, 0);
98 half1Normals[edgei] = half0Normals[edgei];
100 else if (
mag(magLe - nbrMagLe)/avLe > matchTol_)
109 half0Normals[edgei] /= magLe;
110 half1Normals[edgei] /= nbrMagLe;
115 if (maxMatchError > matchTol_)
117 label nbrEdgei = errorEdge + sizeby2;
118 scalar magLe =
mag(half0Normals[errorEdge]);
119 scalar nbrMagLe =
mag(half1Normals[errorEdge]);
120 scalar avLe = 0.5*(magLe + nbrMagLe);
123 <<
"edge " << errorEdge
124 <<
" area does not match neighbour "
125 << nbrEdgei <<
" by "
126 << 100*
mag(magLe - nbrMagLe)/avLe
127 <<
"% -- possible edge ordering problem." <<
endl
128 <<
"patch:" <<
name()
129 <<
" my area:" << magLe
130 <<
" neighbour area:" << nbrMagLe
131 <<
" matching tolerance:" << matchTol_
133 <<
"Mesh edge:" << start() + errorEdge
135 <<
"Neighbour edge:" << start() + nbrEdgei
137 <<
"Other errors also exist, only the largest is reported. "
138 <<
"Please rerun with cyclic debug flag set"
154 if (forwardT().size() > 1 || reverseT().size() > 1)
157 <<
"Transformation tensor is not constant for the cyclic "
158 <<
"patch. Please reconsider your setup and definition of "
159 <<
"cyclic boundaries." <<
endl;
171 const label sizeby2 = deltas.size()/2;
173 scalar maxMatchError = 0;
174 label errorEdge = -1;
176 for (label edgei = 0; edgei < sizeby2; ++edgei)
178 scalar avL = 0.5*(magL[edgei] + magL[edgei + sizeby2]);
182 mag(magL[edgei] - magL[edgei + sizeby2])/avL
191 mag(magL[edgei] - magL[edgei + sizeby2])/avL
197 scalar di = deltas[edgei];
198 scalar dni = deltas[edgei + sizeby2];
200 w[edgei] = dni/(di + dni);
201 w[edgei + sizeby2] = 1 - w[edgei];
205 if (maxMatchError > matchTol_)
207 scalar avL = 0.5*(magL[errorEdge] + magL[errorEdge + sizeby2]);
210 <<
"edge " << errorEdge <<
" and " << errorEdge + sizeby2
211 <<
" areas do not match by "
212 << 100*
mag(magL[errorEdge] - magL[errorEdge + sizeby2])/avL
213 <<
"% -- possible edge ordering problem." <<
nl
214 <<
"Cyclic area match tolerance = "
231 const label sizeby2 = deltas.size()/2;
233 for (label edgei = 0; edgei < sizeby2; ++edgei)
235 scalar di = deltas[edgei];
236 scalar dni = deltas[edgei + sizeby2];
238 dc[edgei] = 1.0/(di + dni);
239 dc[edgei + sizeby2] = dc[edgei];
281 const label sizeby2 = patchD.size()/2;
284 auto& pdv = tpdv.ref();
289 for (label edgei = 0; edgei < sizeby2; ++edgei)
291 const vector& ddi = patchD[edgei];
292 const vector& dni = patchD[edgei + sizeby2];
294 pdv[edgei] = ddi - dni;
295 pdv[edgei + sizeby2] = -pdv[edgei];
300 for (label edgei = 0; edgei < sizeby2; ++edgei)
302 const vector& ddi = patchD[edgei];
303 const vector& dni = patchD[edgei + sizeby2];
305 pdv[edgei] = ddi -
transform(forwardT()[0], dni);
306 pdv[edgei + sizeby2] = -
transform(reverseT()[0], pdv[edgei]);
342 auto& pnf = tpnf.ref();
344 const label sizeby2 = this->size()/2;
346 for (label edgei=0; edgei < sizeby2; ++edgei)
348 pnf[edgei] = interfaceData[edgei + sizeby2];
349 pnf[edgei + sizeby2] = interfaceData[edgei];
374 auto& pnf = tpnf.ref();
376 const label sizeby2 = this->size()/2;
378 for (label edgei=0; edgei < sizeby2; ++edgei)
380 pnf[edgei] = iF[edgeCells[edgei + sizeby2]];
381 pnf[edgei + sizeby2] = iF[edgeCells[edgei]];
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Info<< " "<< writer.output().name()<< nl;}{ const Field< vector > edgeCentres(faMeshTools::flattenEdgeField(aMesh.edgeCentres(), true))
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void size(const label n)
Older name for setAddressableSize.
commsTypes
Communications types.
coupledFaPatch is an abstract base class for patches that couple regions of the computational domain ...
bool parallel() const
Are the cyclic planes parallel.
virtual const labelUList & faceCells() const
Return faceCell addressing: lduInterface virtual function.
coupledFaPatch(const word &name, const labelUList &edgeLabels, const label index, const faBoundaryMesh &bm, const label nbrPolyPatchIndex, const word &patchType)
Construct from components.
virtual tmp< vectorField > delta() const =0
Return delta (P to N) vectors across coupled patch.
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
void makeWeights(scalarField &) const
Make patch weighting factors.
cyclicFaPatch(const word &name, const dictionary &dict, const label index, const faBoundaryMesh &bm, const word &patchType)
Construct from dictionary.
void makeLPN(scalarField &) const
Make patch geodesic distance between P and N.
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
void makeDeltaCoeffs(scalarField &) const
Make patch face - neighbour cell distances.
virtual tmp< labelField > transfer(const Pstream::commsTypes commsType, const labelUList &interfaceData) const
Transfer and return neighbour field.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
static const scalar matchTol_
Relative tolerance (for geometric matching). Is factor of.
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Return the values of the given internal data adjacent to the interface as a field.
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patches after moving points.
virtual const tensorField & forwardT() const
Return face transformation tensor.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Finite area boundary mesh, which is a faPatch list with registered IO, a reference to the associated ...
Finite area patch class. Used for 2-D non-Euclidian finite area method.
virtual label size() const
Patch size is the number of edge labels, but can be overloaded.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
const scalarField & magEdgeLengths() const
Return edge length magnitudes, like the faMesh::magLe() method.
const scalarField & lPN() const
Return patch geodesic distance between P and N.
const labelUList & edgeFaces() const
Return edge-face addressing.
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patch after moving points.
tmp< vectorField > edgeNormals() const
Return edge unit normals, like the faMesh::unitLe() method.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
void patchInternalField(const UList< Type > &internalData, const labelUList &addressing, UList< Type > &pfld) const
Extract internal field next to patch using specified addressing.
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
label index() const noexcept
The index of this patch in the boundaryMesh.
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
UList< label > labelUList
A UList of labels.
constexpr char nl
The newline '\n' character (0x0a).