coupledPolyPatch.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "coupledPolyPatch.H"
29 #include "ListOps.H"
30 #include "transform.H"
31 #include "OFstream.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(coupledPolyPatch, 0);
38 
39  const scalar coupledPolyPatch::defaultMatchTol_ = 1e-4;
40 }
41 
42 
43 const Foam::Enum
44 <
46 >
48 ({
49  { transformType::UNKNOWN, "unknown" },
50  { transformType::ROTATIONAL, "rotational" },
51  { transformType::TRANSLATIONAL, "translational" },
52  { transformType::COINCIDENTFULLMATCH, "coincidentFullMatch" },
53  { transformType::NOORDERING, "noOrdering" },
54 });
55 
56 
57 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
58 
60 {
61  os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
62 }
63 
64 
66 (
67  Ostream& os,
68  const pointField& points,
69  const labelList& pointLabels
70 )
71 {
73  {
74  writeOBJ(os, points[pointLabels[i]]);
75  }
76 }
77 
78 
80 (
81  Ostream& os,
82  const point& p0,
83  const point& p1,
84  label& vertI
85 )
86 {
87  writeOBJ(os, p0);
88  vertI++;
89 
90  writeOBJ(os, p1);
91  vertI++;
92 
93  os << "l " << vertI-1 << ' ' << vertI << nl;
94 }
95 
96 
98 (
99  const fileName& fName,
100  const UList<face>& faces,
101  const pointField& points
102 )
103 {
104  OFstream os(fName);
105 
106  Map<label> foamToObj(4*faces.size());
107 
108  label vertI = 0;
109 
110  forAll(faces, i)
111  {
112  const face& f = faces[i];
113 
114  forAll(f, fp)
115  {
116  if (foamToObj.insert(f[fp], vertI))
117  {
118  writeOBJ(os, points[f[fp]]);
119  vertI++;
120  }
121  }
122 
123  os << 'l';
124  forAll(f, fp)
125  {
126  os << ' ' << foamToObj[f[fp]]+1;
127  }
128  os << ' ' << foamToObj[f[0]]+1 << nl;
129  }
130 }
131 
132 
134 (
135  const UList<face>& faces,
136  const pointField& points,
138 )
139 {
140  pointField anchors(faces.size());
141 
142  if (transform != COINCIDENTFULLMATCH)
143  {
144  // Return the first point
145  forAll(faces, facei)
146  {
147  anchors[facei] = points[faces[facei][0]];
148  }
149  }
150  else
151  {
152  // Make anchor point unique
153  forAll(faces, facei)
154  {
155  const face& f = faces[facei];
156 
157  bool unique = true;
158 
159  forAll(f, fp1)
160  {
161  const point& p1 = points[f[fp1]];
162 
163  unique = true;
164 
165  for (label fp2 = 0; fp2 < f.size(); ++fp2)
166  {
167  if (f[fp1] == f[fp2])
168  {
169  continue;
170  }
171 
172  const point& p2 = points[f[fp2]];
173 
174  // TODO: Change to a tolerance and possibly select closest
175  // point to the origin
176  if (p1 == p2)
177  {
178  unique = false;
179  break;
180  }
181  }
182 
183  if (unique)
184  {
185  anchors[facei] = p1;
186  break;
187  }
188  }
189 
190  if (!unique)
191  {
192  anchors[facei] = points[faces[facei][0]];
193  }
194  }
195  }
196 
197  return anchors;
198 }
199 
200 
202 (
203  const UList<face>& faces,
204  const pointField& points,
205  const pointField& faceCentres
206 )
207 {
208  // Calculate typical distance per face
209  scalarField tols(faces.size());
210 
211  forAll(faces, facei)
212  {
213  const point& cc = faceCentres[facei];
214 
215  const face& f = faces[facei];
216 
217  // 1. calculate a typical size of the face. Use maximum distance
218  // to face centre
219  scalar maxLenSqr = -GREAT;
220  // 2. as measure of truncation error when comparing two coordinates
221  // use SMALL * maximum component
222  scalar maxCmpt = -GREAT;
223 
224  forAll(f, fp)
225  {
226  const point& pt = points[f[fp]];
227  maxLenSqr = max(maxLenSqr, magSqr(pt - cc));
228  maxCmpt = max(maxCmpt, cmptMax(cmptMag(pt)));
229  }
230 
231  tols[facei] = max
232  (
233  SMALL,
234  max(SMALL*maxCmpt, Foam::sqrt(maxLenSqr))
235  );
236  }
237  return tols;
238 }
239 
240 
242 (
243  const pointField& points,
244  const face& f,
245  const point& anchor,
246  const scalar tol
247 )
248 {
249  label anchorFp = -1;
250  scalar minDistSqr = GREAT;
251 
252  forAll(f, fp)
253  {
254  scalar distSqr = magSqr(anchor - points[f[fp]]);
255 
256  if (distSqr < minDistSqr)
257  {
258  minDistSqr = distSqr;
259  anchorFp = fp;
260  }
261  }
262 
263  if (anchorFp == -1 || Foam::sqrt(minDistSqr) > tol)
264  {
265  return -1;
266  }
267  else
268  {
269  // Check that anchor is unique.
270  forAll(f, fp)
271  {
272  scalar distSqr = magSqr(anchor - points[f[fp]]);
273 
274  if (distSqr == minDistSqr && fp != anchorFp)
275  {
277  << "Cannot determine unique anchor point on face "
279  << endl
280  << "Both at index " << anchorFp << " and " << fp
281  << " the vertices have the same distance "
282  << Foam::sqrt(minDistSqr)
283  << " to the anchor " << anchor
284  << ". Continuing but results might be wrong."
285  << nl << endl;
286  }
287  }
288 
289  // Positive rotation
290  return (f.size() - anchorFp) % f.size();
291  }
292 }
293 
294 
296 (
297  const vectorField& Cf,
298  const vectorField& Cr,
299  const vectorField& nf,
300  const vectorField& nr,
301  const scalarField& smallDist,
302  const scalar absTol,
304 ) const
305 {
306  if (debug)
307  {
308  Pout<< "coupledPolyPatch::calcTransformTensors : " << name() << endl
309  << " transform:" << transformTypeNames[transform] << nl
310  << " (half)size:" << Cf.size() << nl
311  << " absTol:" << absTol << nl
312  << " smallDist min:" << min(smallDist) << nl
313  << " smallDist max:" << max(smallDist) << nl
314  << " sum(mag(nf & nr)):" << sum(mag(nf & nr)) << endl;
315  }
316 
317  // Tolerance calculation.
318  // - normal calculation: assume absTol is the absolute error in a
319  // single normal/transformation calculation. Consists both of numerical
320  // precision (on the order of SMALL and of writing precision
321  // (from e.g. decomposition)
322  // Then the overall error of summing the normals is sqrt(size())*absTol
323  // - separation calculation: pass in from the outside an allowable error.
324 
325  if (Cf.size() == 0)
326  {
327  // Dummy geometry. Assume non-separated, parallel.
328  separation_.setSize(0);
329  forwardT_.clear();
330  reverseT_.clear();
331  collocated_.setSize(0);
332  }
333  else
334  {
335  scalar error = absTol*Foam::sqrt(1.0*Cf.size());
336 
337  if (debug)
338  {
339  Pout<< " error:" << error << endl;
340  }
341 
342  if
343  (
344  transform == ROTATIONAL
345  || (
346  transform != TRANSLATIONAL
347  && transform != COINCIDENTFULLMATCH
348  && (sum(mag(nf & nr)) < Cf.size() - error)
349  )
350  )
351  {
352  // Type is rotation or unknown and normals not aligned
353 
354  // Assume per-face differing transformation, correct later
355 
356  separation_.setSize(0);
357 
358  forwardT_.setSize(Cf.size());
359  reverseT_.setSize(Cf.size());
360  collocated_.setSize(Cf.size());
361  collocated_ = false;
362 
363  forAll(forwardT_, facei)
364  {
365  forwardT_[facei] = rotationTensor(-nr[facei], nf[facei]);
366  reverseT_[facei] = rotationTensor(nf[facei], -nr[facei]);
367  }
368 
369  if (debug)
370  {
371  Pout<< " sum(mag(forwardT_ - forwardT_[0])):"
372  << sum(mag(forwardT_ - forwardT_[0]))
373  << endl;
374  }
375 
376  if (sum(mag(forwardT_ - forwardT_[0])) < error)
377  {
378  forwardT_.setSize(1);
379  reverseT_.setSize(1);
380  collocated_.setSize(1);
381 
382  if (debug)
383  {
384  Pout<< " difference in rotation less than"
385  << " local tolerance "
386  << error << ". Assuming uniform rotation." << endl;
387  }
388  }
389  }
390  else
391  {
392  // Translational or (unknown and normals aligned)
393 
394  forwardT_.setSize(0);
395  reverseT_.setSize(0);
396 
397  separation_ = Cr - Cf;
398 
399  collocated_.setSize(separation_.size());
400 
401  // Three situations:
402  // - separation is zero. No separation.
403  // - separation is same. Single separation vector.
404  // - separation differs per face. Separation vectorField.
405 
406  // Check for different separation per face
407  bool sameSeparation = true;
408  bool doneWarning = false;
409 
410  forAll(separation_, facei)
411  {
412  scalar smallSqr = sqr(smallDist[facei]);
413 
414  collocated_[facei] = (magSqr(separation_[facei]) < smallSqr);
415 
416  // Check if separation differing w.r.t. face 0.
417  if (magSqr(separation_[facei] - separation_[0]) > smallSqr)
418  {
419  sameSeparation = false;
420 
421  if (!doneWarning && debug)
422  {
423  doneWarning = true;
424 
425  Pout<< " separation " << separation_[facei]
426  << " at " << facei
427  << " differs from separation[0] " << separation_[0]
428  << " by more than local tolerance "
429  << smallDist[facei]
430  << ". Assuming non-uniform separation." << endl;
431  }
432  }
433  }
434 
435  if (sameSeparation)
436  {
437  // Check for zero separation (at 0 so everywhere)
438  if (collocated_[0])
439  {
440  if (debug)
441  {
442  Pout<< " separation " << mag(separation_[0])
443  << " less than local tolerance " << smallDist[0]
444  << ". Assuming zero separation." << endl;
445  }
446 
447  separation_.setSize(0);
448  collocated_ = boolList(1, true);
449  }
450  else
451  {
452  if (debug)
453  {
454  Pout<< " separation " << mag(separation_[0])
455  << " more than local tolerance " << smallDist[0]
456  << ". Assuming uniform separation." << endl;
457  }
458 
459  separation_.setSize(1);
460  collocated_ = boolList(1, false);
461  }
462  }
463  }
464  }
465 
466  if (debug)
467  {
468  Pout<< " separation_:" << separation_.size() << nl
469  << " forwardT size:" << forwardT_.size() << endl;
470  }
471 }
472 
473 
474 // * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * * * //
475 
477 (
478  const word& name,
479  const label size,
480  const label start,
481  const label index,
482  const polyBoundaryMesh& bm,
483  const word& patchType,
485 )
486 :
487  polyPatch(name, size, start, index, bm, patchType),
488  matchTolerance_(defaultMatchTol_),
489  transform_(transform)
490 {}
491 
492 
494 (
495  const word& name,
496  const dictionary& dict,
497  const label index,
498  const polyBoundaryMesh& bm,
499  const word& patchType
500 )
501 :
502  polyPatch(name, dict, index, bm, patchType),
503  matchTolerance_(dict.lookupOrDefault("matchTolerance", defaultMatchTol_)),
504  transform_
505  (
506  transformTypeNames.lookupOrDefault
507  (
508  "transform",
509  dict,
510  transformType::UNKNOWN
511  )
512  )
513 {}
514 
515 
517 (
518  const coupledPolyPatch& pp,
519  const polyBoundaryMesh& bm
520 )
521 :
522  polyPatch(pp, bm),
523  matchTolerance_(pp.matchTolerance_),
524  transform_(pp.transform_)
525 {}
526 
527 
529 (
530  const coupledPolyPatch& pp,
531  const polyBoundaryMesh& bm,
532  const label index,
533  const label newSize,
534  const label newStart
535 )
536 :
537  polyPatch(pp, bm, index, newSize, newStart),
538  matchTolerance_(pp.matchTolerance_),
539  transform_(pp.transform_)
540 {}
541 
542 
544 (
545  const coupledPolyPatch& pp,
546  const polyBoundaryMesh& bm,
547  const label index,
548  const labelUList& mapAddressing,
549  const label newStart
550 )
551 :
552  polyPatch(pp, bm, index, mapAddressing, newStart),
553  matchTolerance_(pp.matchTolerance_),
554  transform_(pp.transform_)
555 {}
556 
557 
558 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
559 
561 {}
562 
563 
564 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
565 
567 {
568  polyPatch::write(os);
569  //if (matchTolerance_ != defaultMatchTol_)
570  {
571  os.writeEntry("matchTolerance", matchTolerance_);
572  os.writeEntry("transform", transformTypeNames[transform_]);
573  }
574 }
575 
576 
577 // ************************************************************************* //
Foam::coupledPolyPatch::getAnchorPoints
static pointField getAnchorPoints(const UList< face > &, const pointField &, const transformType)
Get a unique anchor point for all faces.
Definition: coupledPolyPatch.C:134
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:78
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:51
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Foam::polyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: polyPatch.C:403
Foam::coupledPolyPatch::transformTypeNames
static const Enum< transformType > transformTypeNames
Definition: coupledPolyPatch.H:69
Foam::Map< label >
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:72
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::coupledPolyPatch
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
Definition: coupledPolyPatch.H:53
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:90
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
OFstream.H
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:400
Foam::cmptMax
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:253
coupledPolyPatch.H
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< vector >
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::coupledPolyPatch::getRotation
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.
Definition: coupledPolyPatch.C:242
Foam::coupledPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: coupledPolyPatch.C:566
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::coupledPolyPatch::calcFaceTol
static scalarField calcFaceTol(const UList< face > &faces, const pointField &points, const pointField &faceCentres)
Calculate typical tolerance per face. Is currently max distance.
Definition: coupledPolyPatch.C:202
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:99
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:84
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::nl
constexpr char nl
Definition: Ostream.H:372
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::coupledPolyPatch::~coupledPolyPatch
virtual ~coupledPolyPatch()
Destructor.
Definition: coupledPolyPatch.C:560
Foam::List< label >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::start
label ListType::const_reference const label start
Definition: ListOps.H:408
Foam::UList< face >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::coupledPolyPatch::calcTransformTensors
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.
Definition: coupledPolyPatch.C:296
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
Foam::UList::size
void size(const label n) noexcept
Override size to be inconsistent with allocated storage.
Definition: UListI.H:360
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:109
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
ListOps.H
Various functions to operate on Lists.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
p0
const volScalarField & p0
Definition: EEqn.H:36
transform.H
3D tensor transformation operations.
Foam::rotationTensor
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:49
Foam::coupledPolyPatch::writeOBJ
static void writeOBJ(Ostream &os, const point &pt)
Write point in OBJ format.
Definition: coupledPolyPatch.C:59
Foam::coupledPolyPatch::transformType
transformType
Definition: coupledPolyPatch.H:59
Foam::coupledPolyPatch::coupledPolyPatch
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.
Definition: coupledPolyPatch.C:477
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::error
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:64
pointLabels
labelList pointLabels(nPoints, -1)