globalIndexAndTransform.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 
27 #include "cyclicPolyPatch.H"
28 #include "DynamicField.H"
29 #include "globalMeshData.H"
30 
31 // * * * * * * * * * * * * Private Static Data Members * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(globalIndexAndTransform, 0);
36 }
37 
38 
39 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 
41 Foam::label Foam::globalIndexAndTransform::matchTransform
42 (
43  const List<vectorTensorTransform>& refTransforms,
44  label& matchedRefTransformI,
45  const vectorTensorTransform& testTransform,
46  scalar tolerance,
47  bool checkBothSigns
48 ) const
49 {
50  matchedRefTransformI = -1;
51 
52  forAll(refTransforms, i)
53  {
54  const vectorTensorTransform& refTransform = refTransforms[i];
55 
56  scalar maxVectorMag = sqrt
57  (
58  max(magSqr(testTransform.t()), magSqr(refTransform.t()))
59  );
60 
61  // Test the difference between vector parts to see if it is
62  // less than tolerance times the larger vector part magnitude.
63 
64  scalar vectorDiff =
65  mag(refTransform.t() - testTransform.t())
66  /(maxVectorMag + vSmall)
67  /tolerance;
68 
69  // Test the difference between tensor parts to see if it is
70  // less than the tolerance. sqrt(3.0) factor used to scale
71  // differences as this is magnitude of a rotation tensor. If
72  // neither transform has a rotation, then the test is not
73  // necessary.
74 
75  scalar tensorDiff = 0;
76 
77  if (refTransform.hasR() || testTransform.hasR())
78  {
79  tensorDiff =
80  mag(refTransform.R() - testTransform.R())
81  /sqrt(3.0)
82  /tolerance;
83  }
84 
85  // ...Diff result is < 1 if the test part matches the ref part
86  // within tolerance
87 
88  if (vectorDiff < 1 && tensorDiff < 1)
89  {
90  matchedRefTransformI = i;
91 
92  return +1;
93  }
94 
95  if (checkBothSigns)
96  {
97  // Test the inverse transform differences too
98 
99  vectorDiff =
100  mag(refTransform.t() + testTransform.t())
101  /(maxVectorMag + vSmall)
102  /tolerance;
103 
104  tensorDiff = 0;
105 
106  if (refTransform.hasR() || testTransform.hasR())
107  {
108  tensorDiff =
109  mag(refTransform.R() - testTransform.R().T())
110  /sqrt(3.0)
111  /tolerance;
112  }
113 
114  if (vectorDiff < 1 && tensorDiff < 1)
115  {
116  matchedRefTransformI = i;
117 
118  return -1;
119  }
120  }
121  }
122 
123  return 0;
124 }
125 
126 
127 void Foam::globalIndexAndTransform::determineTransforms()
128 {
129  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
130 
131  DynamicList<vectorTensorTransform> localTransforms;
132  DynamicField<scalar> localTols;
133 
134  label dummyMatch = -1;
135 
136  forAll(patches, patchi)
137  {
138  const polyPatch& pp = patches[patchi];
139 
140  // Note: special check for unordered cyclics. These are in fact
141  // transform bcs and should probably be split off.
142  // Note: We don't want to be finding transforms for patches marked as
143  // coincident full match. These should have no transform by definition.
144  if
145  (
146  isA<coupledPolyPatch>(pp)
147  && !(
148  isA<cyclicPolyPatch>(pp)
149  && refCast<const cyclicPolyPatch>(pp).transform()
151  )
152  && !(
153  refCast<const coupledPolyPatch>(pp).transform()
155  )
156  )
157  {
158  const coupledPolyPatch& cpp = refCast<const coupledPolyPatch>(pp);
159 
160  if (cpp.separated())
161  {
162  const vectorField& sepVecs = cpp.separation();
163 
164  forAll(sepVecs, sVI)
165  {
166  const vector& sepVec = sepVecs[sVI];
167 
168  if (mag(sepVec) > small)
169  {
170  vectorTensorTransform transform(sepVec);
171 
172  if
173  (
174  matchTransform
175  (
176  localTransforms,
177  dummyMatch,
178  transform,
179  cpp.matchTolerance(),
180  false
181  ) == 0
182  )
183  {
184  localTransforms.append(transform);
185  localTols.append(cpp.matchTolerance());
186  }
187  }
188  }
189  }
190  else if (!cpp.parallel())
191  {
192  const tensorField& transTensors = cpp.reverseT();
193 
194  forAll(transTensors, tTI)
195  {
196  const tensor& transT = transTensors[tTI];
197 
198  if (mag(transT - I) > small)
199  {
200  vectorTensorTransform transform(transT);
201 
202  if
203  (
204  matchTransform
205  (
206  localTransforms,
207  dummyMatch,
208  transform,
209  cpp.matchTolerance(),
210  false
211  ) == 0
212  )
213  {
214  localTransforms.append(transform);
215  localTols.append(cpp.matchTolerance());
216  }
217  }
218  }
219  }
220  }
221  }
222 
223 
224  // Collect transforms on master
225  List<List<vectorTensorTransform>> allTransforms(Pstream::nProcs());
226  allTransforms[Pstream::myProcNo()] = localTransforms;
227  Pstream::gatherList(allTransforms);
228 
229  // Collect matching tolerance on master
230  List<scalarField> allTols(Pstream::nProcs());
231  allTols[Pstream::myProcNo()] = localTols;
232  Pstream::gatherList(allTols);
233 
234  if (Pstream::master())
235  {
236  localTransforms.clear();
237 
238  forAll(allTransforms, proci)
239  {
240  const List<vectorTensorTransform>& procTransVecs =
241  allTransforms[proci];
242 
243  forAll(procTransVecs, pSVI)
244  {
245  const vectorTensorTransform& transform = procTransVecs[pSVI];
246 
247  if (mag(transform.t()) > small || transform.hasR())
248  {
249  if
250  (
251  matchTransform
252  (
253  localTransforms,
254  dummyMatch,
255  transform,
256  allTols[proci][pSVI],
257  true
258  ) == 0
259  )
260  {
261  localTransforms.append(transform);
262  }
263  }
264  }
265  }
266  }
267 
268  transforms_.transfer(localTransforms);
269  Pstream::scatter(transforms_);
270 }
271 
272 
273 void Foam::globalIndexAndTransform::determineTransformPermutations()
274 {
275  label nTransformPermutations = pow(label(3), transforms_.size());
276 
277  transformPermutations_.setSize(nTransformPermutations);
278 
279  forAll(transformPermutations_, tPI)
280  {
281  vectorTensorTransform transform;
282 
283  label transformIndex = tPI;
284 
285  // Invert the ternary index encoding using repeated division by
286  // three
287 
288  forAll(transforms_, b)
289  {
290  const label w = (transformIndex % 3) - 1;
291 
292  transformIndex /= 3;
293 
294  if (w > 0)
295  {
296  transform &= transforms_[b];
297  }
298  else if (w < 0)
299  {
300  transform &= inv(transforms_[b]);
301  }
302  }
303 
304  transformPermutations_[tPI] = transform;
305  }
306 
307 
308  // Encode index for 0 sign
309  labelList permutationIndices(nIndependentTransforms(), 0);
310  nullTransformIndex_ = encodeTransformIndex(permutationIndices);
311 }
312 
313 
314 void Foam::globalIndexAndTransform::determinePatchTransformSign()
315 {
316  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
317 
318  patchTransformSign_.setSize(patches.size(), labelPair(-1, 0));
319 
320  forAll(patches, patchi)
321  {
322  const polyPatch& pp = patches[patchi];
323 
324  // Note: special check for unordered cyclics. These are in fact
325  // transform bcs and should probably be split off.
326  // Note: We don't want to be finding transforms for patches marked as
327  // coincident full match. These should have no transform by definition.
328  if
329  (
330  isA<coupledPolyPatch>(pp)
331  && !(
332  isA<cyclicPolyPatch>(pp)
333  && refCast<const cyclicPolyPatch>(pp).transform()
335  )
336  && !(
337  refCast<const coupledPolyPatch>(pp).transform()
339  )
340  )
341  {
342  const coupledPolyPatch& cpp = refCast<const coupledPolyPatch>(pp);
343 
344  if (cpp.separated())
345  {
346  const vectorField& sepVecs = cpp.separation();
347 
348  // This loop is implicitly expecting only a single
349  // value for separation()
350  forAll(sepVecs, sVI)
351  {
352  const vector& sepVec = sepVecs[sVI];
353 
354  if (mag(sepVec) > small)
355  {
356  vectorTensorTransform t(sepVec);
357 
358  label matchTransI;
359  label sign = matchTransform
360  (
361  transforms_,
362  matchTransI,
363  t,
364  cpp.matchTolerance(),
365  true
366  );
367  patchTransformSign_[patchi] =
368  labelPair(matchTransI, sign);
369  }
370  }
371 
372  }
373  else if (!cpp.parallel())
374  {
375  const tensorField& transTensors = cpp.reverseT();
376 
377  // This loop is implicitly expecting only a single
378  // value for reverseT()
379  forAll(transTensors, tTI)
380  {
381  const tensor& transT = transTensors[tTI];
382 
383  if (mag(transT - I) > small)
384  {
385  vectorTensorTransform t(transT);
386 
387  label matchTransI;
388  label sign = matchTransform
389  (
390  transforms_,
391  matchTransI,
392  t,
393  cpp.matchTolerance(),
394  true
395  );
396 
397  patchTransformSign_[patchi] =
398  labelPair(matchTransI, sign);
399  }
400  }
401  }
402  }
403  }
404 }
405 
406 
407 bool Foam::globalIndexAndTransform::uniqueTransform
408 (
409  const point& pt,
410  labelPairList& trafos,
411  const label patchi,
412  const labelPair& patchTrafo
413 ) const
414 {
415  if (findIndex(trafos, patchTrafo) == -1)
416  {
417  // New transform. Check if already have 3
418  if (trafos.size() == 3)
419  {
420  if (patchi > -1)
421  {
423  << "Point " << pt
424  << " is on patch " << mesh_.boundaryMesh()[patchi].name();
425  }
426  else
427  {
429  << "Point " << pt << " is on a coupled patch";
430  }
431  Warning
432  << " with transformation " << patchTrafo
433  << " but also on 3 other patches with transforms "
434  << trafos << nl
435  << "This is not a space filling tiling and might"
436  << " indicate a setup problem and give problems"
437  << " for e.g. lagrangian tracking or interpolation" << endl;
438 
439  // Already warned so no need to extend more
440  trafos.clear();
441  return false;
442  }
443 
444  return true;
445  }
446  else
447  {
448  return false;
449  }
450 }
451 
452 
453 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
454 
455 Foam::globalIndexAndTransform::globalIndexAndTransform(const polyMesh& mesh)
456 :
457  mesh_(mesh),
458  transforms_(),
459  transformPermutations_(),
460  patchTransformSign_()
461 {
462  determineTransforms();
463 
464  determineTransformPermutations();
465 
466  determinePatchTransformSign();
467 
468  if (debug && transforms_.size() > 0)
469  {
470  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
471 
472  Info<< "Determined global transforms :" << endl;
473  Info<< "\t\ttranslation\trotation" << endl;
474  forAll(transforms_, i)
475  {
476  Info<< '\t' << i << '\t';
477  const vectorTensorTransform& trafo = transforms_[i];
478  if (trafo.hasR())
479  {
480  Info<< trafo.t() << '\t' << trafo.R();
481  }
482  else
483  {
484  Info<< trafo.t() << '\t' << "---";
485  }
486  Info<< endl;
487  }
488  Info<< endl;
489 
490 
491  Info<< "\tpatch\ttransform\tsign" << endl;
492  forAll(patchTransformSign_, patchi)
493  {
494  if (patchTransformSign_[patchi].first() != -1)
495  {
496  Info<< '\t' << patches[patchi].name()
497  << '\t' << patchTransformSign_[patchi].first()
498  << '\t' << patchTransformSign_[patchi].second()
499  << endl;
500  }
501  }
502  Info<< endl;
503 
504 
505  Info<< "Permutations of transformations:" << endl
506  << "\t\ttranslation\trotation" << endl;
507  forAll(transformPermutations_, i)
508  {
509  Info<< '\t' << i << '\t';
510  const vectorTensorTransform& trafo = transformPermutations_[i];
511  if (trafo.hasR())
512  {
513  Info<< trafo.t() << '\t' << trafo.R();
514  }
515  else
516  {
517  Info<< trafo.t() << '\t' << "---";
518  }
519  Info<< endl;
520  }
521  Info<< "nullTransformIndex:" << nullTransformIndex() << endl
522  << endl;
523  }
524 
525 
526  if (transforms_.size() > 0)
527  {
528  // Check that the transforms are space filling : any point
529  // can only use up to three transforms
530 
531  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
532 
533 
534  // 1. Collect transform&sign per point and do local check
535 
536  List<labelPairList> pointToTrafos(mesh_.nPoints());
537 
538  forAll(patches, patchi)
539  {
540  const polyPatch& pp = patches[patchi];
541 
542  const labelPair& transSign = patchTransformSign_[patchi];
543 
544  if (transSign.first() > -1)
545  {
546  const labelList& mp = pp.meshPoints();
547  forAll(mp, i)
548  {
549  labelPairList& trafos = pointToTrafos[mp[i]];
550 
551  bool newTransform = uniqueTransform
552  (
553  mesh_.points()[mp[i]],
554  trafos,
555  patchi,
556  transSign
557  );
558 
559  if (newTransform)
560  {
561  trafos.append(transSign);
562  }
563  }
564  }
565  }
566 
567 
568  // Synchronise across collocated (= untransformed) points
569  // TBD: there is a big problem in that globalMeshData uses
570  // globalIndexAndTransform. Triggers recursion.
571  if (false)
572  {
573  const globalMeshData& gmd = mesh_.globalData();
574  const indirectPrimitivePatch& cpp = gmd.coupledPatch();
575  const labelList& meshPoints = cpp.meshPoints();
576  const mapDistribute& slavesMap = gmd.globalCoPointSlavesMap();
577  const labelListList& slaves = gmd.globalCoPointSlaves();
578 
579  List<labelPairList> elems(slavesMap.constructSize());
580  forAll(meshPoints, i)
581  {
582  elems[i] = pointToTrafos[meshPoints[i]];
583  }
584 
585  // Pull slave data onto master. No need to update transformed slots.
586  slavesMap.distribute(elems, false);
587 
588  // Combine master data with slave data
589  forAll(slaves, i)
590  {
591  labelPairList& trafos = elems[i];
592 
593  const labelList& slavePoints = slaves[i];
594 
595  // Combine master with untransformed slave data
596  forAll(slavePoints, j)
597  {
598  const labelPairList& slaveTrafos = elems[slavePoints[j]];
599 
600  forAll(slaveTrafos, slaveI)
601  {
602  bool newTransform = uniqueTransform
603  (
604  mesh_.points()[meshPoints[i]],
605  trafos,
606  -1,
607  slaveTrafos[slaveI]
608  );
609 
610  if (newTransform)
611  {
612  trafos.append(slaveTrafos[slaveI]);
613  }
614  }
615  }
616  }
617  }
618  }
619 }
620 
621 
622 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
dimensionedScalar sign(const dimensionedScalar &ds)
const mapDistribute & globalCoPointSlavesMap() const
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
const word & name() const
Return name.
Definition: IOobject.H:297
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Vector-tensor class used to perform translations and rotations in 3D space.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:427
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:421
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1003
A list of faces which address into the list of points.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
dynamicFvMesh & mesh
static const Identity< scalar > I
Definition: Identity.H:93
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
label nullTransformIndex() const
Return the transformIndex (index in transformPermutations)
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1174
List< label > labelList
A List of labels.
Definition: labelList.H:56
List< labelPair > labelPairList
List of labelPairs.
Definition: labelPair.H:49
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Foam::polyBoundaryMesh.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:265
defineTypeNameAndDebug(combustionModel, 0)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
messageStream Warning
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label constructSize() const
Constructed data size.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:409
label patchi
Class containing processor-to-processor mapping information.
vector point
Point is a vector.
Definition: point.H:41
#define WarningInFunction
Report a warning using Foam::Warning.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
label nPoints() const
Field< vector > vectorField
Specialisation of Field<T> for vector.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
const Type & first() const
Return first.
Definition: Pair.H:87
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
const dimensionedScalar mp
Proton mass.
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
const labelListList & globalCoPointSlaves() const