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