LagrangianMesh.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) 2025 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 "IOmanip.H"
27 #include "indexedOctree.H"
28 #include "labelIOField.H"
29 #include "LagrangianMesh.H"
30 #include "LagrangianMeshLocation.H"
31 #include "LagrangianModels.H"
32 #include "ListOps.H"
33 #include "meshObjects.H"
34 #include "Time.H"
35 #include "tracking.H"
36 #include "treeDataCell.H"
37 #include "debug.H"
38 
43 
47 
48 #include "polyTopoChangeMap.H"
49 #include "polyMeshMap.H"
50 #include "polyDistributionMap.H"
51 
52 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
57 
58  const word LagrangianMesh::prefix("Lagrangian");
59 
60  const word LagrangianMesh::coordinatesName("coordinates");
61  const word LagrangianMesh::positionName("position");
62  const word LagrangianMesh::stateName("state");
63  const word LagrangianMesh::fractionName("fraction");
64 
67  {"copy", "inPlace"};
68 
72  (
73  (LagrangianMesh::typeName + "Permutation").c_str(),
76  );
77 
78  const NamedEnum<LagrangianMesh::partitioningAlgorithm, 2>
80  {"bin", "quick"};
81 
85  (
86  (LagrangianMesh::typeName + "Partitioning").c_str(),
89  );
90 }
91 
92 
93 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
94 
95 void Foam::LagrangianMesh::printGroups(const bool header) const
96 {
97  checkPtr(offsetsPtr_, "Offsets");
98  const labelList& offsets = offsetsPtr_();
99 
100  // Shorthand casts of group indices
101  static const label completei =
102  static_cast<label>(LagrangianGroup::complete);
103  static const label inInternalMeshi =
105  static const label onPatchZeroi =
106  static_cast<label>(LagrangianGroup::onPatchZero);
107 
108  // Determine the number of non-processor patches
109  label nNonProcPatches = 0;
110  for (; nNonProcPatches < boundary().size(); ++ nNonProcPatches)
111  {
112  const LagrangianPatch& p = boundary()[nNonProcPatches];
113 
114  if (isA<processorLagrangianPatch>(p)) break;
115  }
116 
117  // Determine if there are processor patches and/or if any patches are
118  // referred to by processor-cyclics
119  bool hasProcPatches = false;
120  boolList patchIsReferred(nNonProcPatches, false);
121  for (label patchi = nNonProcPatches; patchi < boundary().size(); ++ patchi)
122  {
123  const LagrangianPatch& p = boundary()[patchi];
124 
125  if (isA<processorCyclicLagrangianPatch>(p))
126  {
127  patchIsReferred
128  [
129  refCast<const processorCyclicLagrangianPatch>(p)
130  .referPatchIndex()
131  ] = true;
132  }
133  else // isA<processorLagrangianPatch>(p)
134  {
135  hasProcPatches = true;
136  }
137  }
138  reduce(hasProcPatches, orOp<bool>());
139  Pstream::listCombineGather(patchIsReferred, orEqOp<bool>());
140  Pstream::listCombineScatter(patchIsReferred);
141 
142  // Get the indices of patches that are referred to by processor-cyclics
143  const labelList referredPatches(findIndices(patchIsReferred, true));
144  const labelList patchReferred(invert(nNonProcPatches, referredPatches));
145 
146  // Determine the indices of the groups in the table after those relating to
147  // the global patches
148  const label onProcPatchi = onPatchZeroi + nNonProcPatches;
149  const label onProcCyclicPatchZeroi = onProcPatchi + hasProcPatches;
150  const label toBeRemovedi = onProcCyclicPatchZeroi + referredPatches.size();
151 
152  // Build the column names
153  wordList columnNames(toBeRemovedi + 1);
154  columnNames[completei] = Foam::name(LagrangianGroup::complete);
155  columnNames[inInternalMeshi] = Foam::name(LagrangianGroup::inInternalMesh);
156  for (label patchi = 0; patchi < nNonProcPatches; ++ patchi)
157  {
158  columnNames[onPatchZeroi + patchi] =
159  isA<nonConformalErrorLagrangianPatch>(boundary()[patchi])
160  || isA<internalLagrangianPatch>(boundary()[patchi])
161  ? word::null
162  : boundary()[patchi].name();
163  }
164  if (hasProcPatches)
165  {
166  columnNames[onProcPatchi] = "(processor)";
167  }
168  forAll(referredPatches, referredPatchi)
169  {
170  columnNames[onProcCyclicPatchZeroi + referredPatchi] =
171  "(processor)"
172  + boundary()[referredPatches[referredPatchi]].name();
173  }
174  columnNames[toBeRemovedi] = Foam::name(LagrangianGroup::toBeRemoved);
175 
176  // Print the column names if the header is requested, and finish
177  if (header)
178  {
179  forAll(columnNames, i)
180  {
181  if (i && columnNames[i].size()) Info<< ' ';
182 
183  Info<< columnNames[i];
184  }
185 
186  return;
187  }
188 
189  // Construct the numbers of elements in each column
190  labelList columnNumbers(toBeRemovedi + 1, 0);
191  columnNumbers[completei] = offsets[completei + 1] - offsets[completei];
192  columnNumbers[inInternalMeshi] =
193  offsets[inInternalMeshi + 1] - offsets[inInternalMeshi];
194  for (label patchi = 0; patchi < nNonProcPatches; ++ patchi)
195  {
196  columnNumbers[onPatchZeroi + patchi] =
197  offsets[onPatchZeroi + patchi + 1] - offsets[onPatchZeroi + patchi];
198  }
199  for (label patchi = nNonProcPatches; patchi < boundary().size(); ++ patchi)
200  {
201  const LagrangianPatch& p = boundary()[patchi];
202 
203  if (isA<processorCyclicLagrangianPatch>(p))
204  {
205  columnNumbers
206  [
207  patchReferred
208  [
209  refCast<const processorCyclicLagrangianPatch>(p)
210  .referPatchIndex()
211  ]
212  ] +=
213  offsets[onPatchZeroi + patchi + 1]
214  - offsets[onPatchZeroi + patchi];
215  }
216  else // isA<processorLagrangianPatch>(p)
217  {
218  columnNumbers[onProcPatchi] +=
219  offsets[onPatchZeroi + patchi + 1]
220  - offsets[onPatchZeroi + patchi];
221  }
222  }
223  Pstream::listCombineGather(columnNumbers, plusEqOp<label>());
224  Pstream::listCombineScatter(columnNumbers);
225 
226  // Print the numbers
227  forAll(columnNames, i)
228  {
229  if (i && columnNames[i].size()) Info<< ' ';
230 
231  const string::size_type columnNumberDigits =
232  Foam::name(columnNumbers[i]).size();
233 
234  if (columnNumberDigits > columnNames[i].size())
235  {
236  Info << string(columnNames[i].size(), '#').c_str();
237  }
238  else
239  {
240  Info << setw(columnNames[i].size()) << columnNumbers[i];
241  }
242  }
243 }
244 
245 
246 Foam::labelList Foam::LagrangianMesh::partitionBin
247 (
248  labelList& offsets,
249  const List<LagrangianState>& states
250 ) const
251 {
252  // Set the zero index so we only partition from the first previously
253  // incomplete element onwards
254  const label i0 = offsets[1];
255 
256  // Sum the numbers of elements in each group and store in the offsets array
257  offsets = 0;
258  for (label i = i0; i < states.size(); ++ i)
259  {
260  const label groupi = stateToGroupi(states[i]);
261  offsets[groupi + 1] ++;
262  }
263 
264  // Cumulative sum, starting at i0, to generate the offsets for the
265  // non-complete subsets of elements
266  offsets[0] = i0;
267  for (label groupi = 0; groupi < nGroups(); ++ groupi)
268  {
269  offsets[groupi + 1] += offsets[groupi];
270  }
271 
272  // Insert each element into the permutation. Increment the offsets to keep
273  // track of the current insertion position within each group.
274  labelList permutation(states.size() - i0);
275  for (label i = i0; i < states.size(); ++ i)
276  {
277  const label groupi = stateToGroupi(states[i]);
278  permutation[offsets[groupi] - i0] = i;
279  ++ offsets[groupi];
280  }
281 
282  // The offsets are now shifted by one. Copy them back a place and set the
283  // first offset to zero, rather than i0, so that the offsets refer to the
284  // entire mesh, not just the non-complete subsets.
285  for (label groupi = nGroups(); groupi > 0; -- groupi)
286  {
287  offsets[groupi] = offsets[groupi - 1];
288  }
289  offsets[0] = 0;
290 
291  // Return the permutation
292  return permutation;
293 }
294 
295 
296 Foam::labelList Foam::LagrangianMesh::partitionQuick
297 (
298  labelList& offsets,
299  const List<LagrangianState>& states
300 ) const
301 {
302  // Set the zero index so we only partition from the first previously
303  // incomplete element onwards
304  const label i0 = offsets[1];
305 
306  // Initialise the offsets, excluding all previously complete elements
307  offsets = -1;
308  offsets[0] = i0;
309  offsets[nGroups()] = states.size();
310 
311  // Initialise a permutation
312  labelList permutation(states.size() - i0);
313  forAll(permutation, i)
314  {
315  permutation[i] = i + i0;
316  }
317 
318  // Compute the number of iterations needed to order the groups
319  label nIterations = 0;
320  for (label groupi = 1; groupi < nGroups(); groupi *= 2)
321  {
322  ++ nIterations;
323  }
324 
325  // Partition
326  label nDivisions = 1;
327  for (label iterationi = 0; iterationi < nIterations; ++ iterationi)
328  {
329  for (label divisioni = 0; divisioni < nDivisions; ++ divisioni)
330  {
331  const label pivot = (2*divisioni + 1)*nGroups()/(2*nDivisions);
332  if (offsets[pivot] != -1)
333  {
334  continue;
335  }
336 
337  const label pivot0 = divisioni*nGroups()/nDivisions;
338  const label pivot1 = (divisioni + 1)*nGroups()/nDivisions;
339  label i = offsets[pivot0] - i0, j = offsets[pivot1] - i0;
340  while (i < j)
341  {
342  while
343  (
344  i < j
345  && stateToGroupi(states[permutation[i]]) < pivot
346  )
347  {
348  ++ i;
349  }
350  while
351  (
352  i < j
353  && stateToGroupi(states[permutation[j - 1]]) >= pivot
354  )
355  {
356  -- j;
357  }
358  if (i < j)
359  {
360  Swap(permutation[i], permutation[j - 1]);
361  }
362  }
363 
364  offsets[pivot] = i + i0;
365  }
366 
367  nDivisions *= 2;
368  }
369 
370  // Re-include previously complete elements into the complete group
371  offsets[0] = 0;
372 
373  // Return the permutation
374  return permutation;
375 }
376 
377 
378 void Foam::LagrangianMesh::permuteAndResizeFields(const labelList& permutation)
379 {
380  wordHashSet permutedFieldNames;
381  #define PERMUTE_TYPE_FIELDS(Type, GeoField) \
382  { \
383  HashTable<GeoField<Type>*> fields \
384  ( \
385  lookupCurrentFields<GeoField<Type>>() \
386  ); \
387  \
388  forAllIter(typename HashTable<GeoField<Type>*>, fields, iter) \
389  { \
390  if (permutedFieldNames.found(iter()->name())) continue; \
391  \
392  permutedFieldNames.insert(iter()->name()); \
393  \
394  permuteList(permutation, iter()->primitiveFieldRef()); \
395  \
396  resizeContainer(iter()->primitiveFieldRef()); \
397  } \
398  }
407  #undef PERMUTE_TYPE_FIELDS
408 }
409 
410 
411 template<class Type>
412 void Foam::LagrangianMesh::permuteList
413 (
414  const labelList& permutation,
415  UList<Type>& list
416 )
417 {
418  switch (permutationAlgorithm_)
419  {
420  case permutationAlgorithm::copy:
421  permuteListCopy(permutation, list);
422  break;
423 
424  case permutationAlgorithm::inPlace:
425  permuteListInPlace(permutation, list);
426  break;
427  }
428 }
429 
430 
431 template<class Type>
432 void Foam::LagrangianMesh::permuteListCopy
433 (
434  const labelList& permutation,
435  UList<Type>& list
436 )
437 {
438  if (permutation.empty()) return;
439 
440  const label i0 = list.size() - permutation.size();
441 
442  SubList<Type>(list, permutation.size(), i0) =
443  List<Type>(UIndirectList<Type>(list, permutation));
444 }
445 
446 
447 template<class Type>
448 void Foam::LagrangianMesh::permuteListInPlace
449 (
450  const labelList& permutation,
451  UList<Type>& list
452 )
453 {
454  if (permutation.empty()) return;
455 
456  const label i0 = list.size() - permutation.size();
457 
458  labelList& permutationRef = const_cast<labelList&>(permutation);
459 
460  label i = 0;
461  Type t = list[i + i0];
462 
463  forAll(permutation, permutationi)
464  {
465  bool end = permutationRef[i] < 0;
466 
467  while (permutationRef[i] < 0) ++ i;
468 
469  if (end) t = list[i + i0];
470 
472  //Swap(list[permutationRef[i]], t);
473 
474  // Reverse
475  list[i + i0] =
476  permutationRef[permutationRef[i] - i0] < 0
477  ? t : list[permutationRef[i]];
478 
479  permutationRef[i] = - permutationRef[i] - 1;
480  i = - permutationRef[i] - 1 - i0;
481  }
482 
483  forAll(permutationRef, permutationi)
484  {
485  permutationRef[permutationi] = - permutationRef[permutationi] - 1;
486  }
487 }
488 
489 
490 template<class Container>
491 void Foam::LagrangianMesh::resizeContainer(Container& container) const
492 {
493  container.resize(offsetsPtr_->last());
494 }
495 
496 
497 Foam::LagrangianSubMesh Foam::LagrangianMesh::append
498 (
500  const labelField& celli,
501  const labelField& facei,
502  const labelField& faceTrii
503 )
504 {
505  clearPosition();
506 
507  const LagrangianSubMesh appendMesh
508  (
509  *this,
511  coordinates.size(),
512  size()
513  );
514 
515  if (statesPtr_.valid())
516  {
517  states().resize(appendMesh.end(), LagrangianState::none);
518  }
519 
520  if (receivePatchFacePtr_.valid())
521  {
522  receivePatchFacePtr_().resize(appendMesh.end(), label(-1));
523  }
524 
525  if (receivePositionPtr_.valid())
526  {
527  receivePositionPtr_().resize(appendMesh.end(), point::nan);
528  }
529 
530  coordinates_.append(coordinates);
531  celli_.append(celli);
532  facei_.append(facei);
533  faceTrii_.append(faceTrii);
534 
535  subAll_.size_ = size();
536 
537  return appendMesh;
538 }
539 
540 
541 Foam::LagrangianSubMesh Foam::LagrangianMesh::append
542 (
543  const labelList& parents
544 )
545 {
546  clearPosition();
547 
548  const LagrangianSubMesh appendMesh
549  (
550  *this,
552  parents.size(),
553  size()
554  );
555 
556  if (statesPtr_.valid())
557  {
558  states().resize(appendMesh.end());
559  appendMesh.sub(static_cast<List<LagrangianState>&>(states())) =
560  UIndirectList<LagrangianState>(states(), parents)();
561  }
562 
563  if (receivePatchFacePtr_.valid())
564  {
565  receivePatchFacePtr_().resize(appendMesh.end());
566  appendMesh.sub(static_cast<List<label>&>(receivePatchFacePtr_())) =
567  UIndirectList<label>(receivePatchFacePtr_(), parents)();
568  }
569 
570  if (receivePositionPtr_.valid())
571  {
572  receivePositionPtr_().resize(appendMesh.end());
573  appendMesh.sub(static_cast<List<point>&>(receivePositionPtr_())) =
574  UIndirectList<point>(receivePositionPtr_(), parents)();
575  }
576 
577  coordinates_.resize(appendMesh.end());
578  appendMesh.sub(static_cast<List<barycentric>&>(coordinates_)) =
579  UIndirectList<barycentric>(coordinates_, parents)();
580  celli_.resize(appendMesh.end());
581  appendMesh.sub(static_cast<labelList&>(celli_)) =
582  UIndirectList<label>(celli_, parents)();
583  facei_.resize(appendMesh.end());
584  appendMesh.sub(static_cast<labelList&>(facei_)) =
585  UIndirectList<label>(facei_, parents)();
586  faceTrii_.resize(appendMesh.end());
587  appendMesh.sub(static_cast<labelList&>(faceTrii_)) =
588  UIndirectList<label>(faceTrii_, parents)();
589 
590  subAll_.size_ = size();
591 
592  return appendMesh;
593 }
594 
595 
596 Foam::wordHashSet Foam::LagrangianMesh::appendSpecifiedFields
597 (
598  const LagrangianSubMesh&
599 ) const
600 {
601  return wordHashSet();
602 }
603 
604 
605 void Foam::LagrangianMesh::injectUnspecifiedFields
606 (
607  const LagrangianInjection& injection,
608  const LagrangianSubMesh& injectionMesh,
609  const wordHashSet& specifiedFieldNames
610 )
611 {
612  // Inject values for unspecified fields using their source conditions
613  wordHashSet injectedFieldNames;
614  #define INJECT_TYPE_FIELDS(Type, GeoField) \
615  { \
616  HashTable<GeoField<Type>*> fields \
617  ( \
618  lookupCurrentFields<GeoField<Type>>() \
619  ); \
620  \
621  forAllIter(typename HashTable<GeoField<Type>*>, fields, iter) \
622  { \
623  if (specifiedFieldNames.found(iter()->name())) continue; \
624  \
625  injectedFieldNames.insert(iter()->name()); \
626  \
627  /* Resize the field */ \
628  iter()->resize(injectionMesh.end()); \
629  \
630  /* Use the source condition to set the new values */ \
631  injectionMesh.sub(*iter()).ref() = \
632  iter()->sources()[injection.name()].value \
633  ( \
634  injection, \
635  injectionMesh \
636  ); \
637  } \
638  }
643  #undef INJECT_TYPE_FIELDS
644 
645  // Special handling for a retained state field
646  #define INJECT_STATE_FIELD(GeoField) \
647  { \
648  if (foundObject<GeoField<label>>(stateName)) \
649  { \
650  GeoField<label>& state = \
651  lookupObjectRef<GeoField<label>>(stateName); \
652  \
653  injectedFieldNames.insert(stateName); \
654  \
655  /* Resize the field */ \
656  state.resize(injectionMesh.end()); \
657  \
658  /* Set unknown state */ \
659  injectionMesh.sub(state).ref() = \
660  static_cast<label>(LagrangianState::none); \
661  } \
662  }
665  #undef INJECT_STATE_FIELD
666 
667  // Make a table of internal fields for which values could not be set
668  wordHashSet internalFieldNames;
669  #define INSERT_INTERNAL_FIELD_NAMES(Type, GeoField) \
670  { \
671  HashTable<GeoField<Type>*> fields \
672  ( \
673  lookupCurrentFields<GeoField<Type>>() \
674  ); \
675  \
676  forAllIter(typename HashTable<GeoField<Type>*>, fields, iter) \
677  { \
678  if (specifiedFieldNames.found(iter()->name())) continue; \
679  if (injectedFieldNames.found(iter()->name())) continue; \
680  \
681  internalFieldNames.insert(iter()->name()); \
682  } \
683  }
688  (
691  );
692  #undef INSERT_INTERNAL_FIELD_NAMES
693 
694  if (!internalFieldNames.empty())
695  {
697  << "Internal fields " << internalFieldNames.sortedToc() << " are "
698  << "registered during an injection event. These fields do not "
699  << "contain source conditions so their new values cannot be "
700  << "assigned. Either specify these fields' values in the calling "
701  << "code, or ensure that they are not registered."
702  << exit(FatalError);
703  }
704 }
705 
706 
707 void Foam::LagrangianMesh::injectUnspecifiedFields
708 (
709  const LagrangianSubMesh& injectionMesh,
710  const wordHashSet& specifiedFieldNames
711 )
712 {
714  #define INSERT_FIELD_NAMES(Type, GeoField) \
715  fieldNames.insert(lookupCurrentFields<GeoField<Type>>().toc());
724  #undef INSERT_FIELD_NAMES
725 
726  if (!fieldNames.empty())
727  {
729  << "Fields " << fieldNames.sortedToc() << " are registered during "
730  << "an injection event. These fields' new values cannot be "
731  << "assigned. Either specify these fields' values in the calling "
732  << "code, or ensure that they are not registered."
733  << exit(FatalError);
734  }
735 }
736 
737 
738 void Foam::LagrangianMesh::birthUnspecifiedFields
739 (
740  const labelList& parents,
741  const LagrangianSubMesh& birthMesh,
742  const wordHashSet& specifiedFieldNames
743 )
744 {
745  wordHashSet birthedFieldNames;
746  #define BIRTH_TYPE_FIELDS(Type, GeoField) \
747  { \
748  HashTable<GeoField<Type>*> fields \
749  ( \
750  lookupCurrentFields<GeoField<Type>>() \
751  ); \
752  \
753  forAllIter(typename HashTable<GeoField<Type>*>, fields, iter) \
754  { \
755  if (specifiedFieldNames.found(iter()->name())) continue; \
756  if (birthedFieldNames.found(iter()->name())) continue; \
757  \
758  birthedFieldNames.insert(iter()->name()); \
759  \
760  /* Resize the field */ \
761  iter()->resize(birthMesh.end()); \
762  \
763  /* Map values from the parent elements */ \
764  birthMesh.sub(*iter()).ref().primitiveFieldRef() = \
765  Field<Type>(UIndirectList<Type>(*iter(), parents)()); \
766  } \
767  }
776  #undef BIRTH_TYPE_FIELDS
777 }
778 
779 
780 void Foam::LagrangianMesh::changer::constructNonConformal() const
781 {
782  bool haveNccPatches = false;
783 
784  forAll(mesh_.boundary(), patchi)
785  {
786  const polyPatch& pp = mesh_.mesh().boundaryMesh()[patchi];
787 
788  if (isA<nonConformalCyclicPolyPatch>(pp))
789  {
790  haveNccPatches = true;
791  break;
792  }
793  }
794 
795  if (!haveNccPatches) return;
796 
797  mesh_.origPatchNccPatchisPtr_.set
798  (
799  new labelListList
800  (
801  mesh_.boundary().size(),
802  labelList()
803  )
804  );
805  mesh_.origPatchNccPatchesPtr_.set
806  (
807  new List<UPtrList<const nonConformalCyclicPolyPatch>>
808  (
809  mesh_.boundary().size(),
810  UPtrList<const nonConformalCyclicPolyPatch>()
811  )
812  );
813  mesh_.nccPatchProcNccPatchisPtr_.set
814  (
815  new labelListList
816  (
817  mesh_.boundary().size(),
819  )
820  );
821 
822  forAll(mesh_.boundary(), patchi)
823  {
824  const polyPatch& pp = mesh_.mesh().boundaryMesh()[patchi];
825 
826  if (isA<nonConformalCyclicPolyPatch>(pp))
827  {
828  const nonConformalCyclicPolyPatch& nccPp =
829  refCast<const nonConformalCyclicPolyPatch>(pp);
830 
831  mesh_.origPatchNccPatchisPtr_()
832  [nccPp.origPatchIndex()].append(patchi);
833  mesh_.origPatchNccPatchesPtr_()
834  [nccPp.origPatchIndex()].append(&nccPp);
835 
836  mesh_.nccPatchProcNccPatchisPtr_()
837  [nccPp.index()][Pstream::myProcNo()] =
838  nccPp.index();
839 
840  if (nccPp.owner()) nccPp.rays();
841  }
842  else if (isA<nonConformalProcessorCyclicPolyPatch>(pp))
843  {
844  const nonConformalProcessorCyclicPolyPatch& ncpcPp =
845  refCast<const nonConformalProcessorCyclicPolyPatch>(pp);
846 
847  mesh_.nccPatchProcNccPatchisPtr_()
848  [ncpcPp.referPatchIndex()][ncpcPp.neighbProcNo()] =
849  ncpcPp.index();
850  }
851  }
852 
853  mesh_.receivePatchFacePtr_.set
854  (
855  new DynamicList<label>(mesh_.size(), label(-1))
856  );
857 
858  mesh_.receivePositionPtr_.set
859  (
860  new DynamicList<point>(mesh_.size(), point::nan)
861  );
862 }
863 
864 
865 void Foam::LagrangianMesh::changer::constructBehind() const
866 {
867  mesh_.fractionBehindPtr_.set
868  (
869  new LagrangianDynamicField<scalar>
870  (
871  IOobject
872  (
873  fractionName + "Behind",
874  mesh_.time().name(),
875  mesh_,
878  ),
879  mesh_,
880  dimensioned<scalar>(dimless, scalar(0)),
881  wordList
882  (
883  mesh_.boundary().size(),
884  calculatedLagrangianPatchScalarField::typeName
885  ),
886  wordList::null(),
887  LagrangianModels::New(mesh_).modelTypeFieldSourceTypes
888  <
889  LagrangianInjection,
890  zeroLagrangianScalarFieldSource,
891  LagrangianSource,
892  internalLagrangianScalarFieldSource
893  >()
894  )
895  );
896 
897  mesh_.nTracksBehindPtr_.set
898  (
899  new LagrangianDynamicField<label>
900  (
901  IOobject
902  (
903  "nTracksBehind",
904  mesh_.time().name(),
905  mesh_,
908  ),
909  mesh_,
910  dimensioned<label>(dimless, label(0)),
911  wordList
912  (
913  mesh_.boundary().size(),
914  calculatedLagrangianPatchLabelField::typeName
915  ),
916  wordList::null(),
917  LagrangianModels::New(mesh_).modelTypeFieldSourceTypes
918  <
919  LagrangianInjection,
920  zeroLagrangianLabelFieldSource,
921  LagrangianSource,
922  internalLagrangianLabelFieldSource
923  >()
924  )
925  );
926 }
927 
928 
929 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
930 
932 (
933  const polyMesh& mesh,
934  const word& name,
937 )
938 :
940  (
941  mesh,
942  name,
943  mesh.boundaryMesh().types(),
944  readOption,
946  )
947 {}
948 
949 
951 (
952  const polyMesh& mesh,
953  const word& name,
954  const wordList& wantedPatchTypes,
957 )
958 :
960  (
961  IOobject
962  (
963  name,
964  mesh.time().name(),
965  prefix,
966  mesh,
967  IOobject::NO_READ,
968  IOobject::AUTO_WRITE
969  )
970  ),
972  mesh_(mesh),
973  coordinates_
974  (
975  IOobject
976  (
977  coordinatesName,
978  time().name(),
979  *this,
980  readOption,
981  IOobject::AUTO_WRITE
982  )
983  ),
984  celli_
985  (
986  IOobject
987  (
988  "cell",
989  time().name(),
990  *this,
991  readOption,
992  IOobject::AUTO_WRITE
993  )
994  ),
995  facei_
996  (
997  IOobject
998  (
999  "face",
1000  time().name(),
1001  *this,
1002  readOption,
1003  IOobject::AUTO_WRITE
1004  )
1005  ),
1006  faceTrii_
1007  (
1008  IOobject
1009  (
1010  "faceTri",
1011  time().name(),
1012  *this,
1013  readOption,
1014  IOobject::AUTO_WRITE
1015  )
1016  ),
1017  boundary_(*this, mesh.boundaryMesh(), wantedPatchTypes),
1018  subAll_
1019  (
1021  (
1022  *this,
1024  size(),
1025  label(0),
1026  label(0)
1027  )
1028  ),
1029  statesPtr_(nullptr),
1030  offsetsPtr_(nullptr),
1031  subMeshIndex_(0),
1032  schemesPtr_(nullptr)
1033 {
1034  writeOpt() = writeOption;
1035 
1036  checkFieldSize(coordinates_);
1037  checkFieldSize(celli_);
1038 
1039  // Read cellFace file if face does not exist. This is useful for testing,
1040  // as it is much easier to write a cellFace file by hand.
1041  if (!celli_.empty() && facei_.empty())
1042  {
1043  labelIOField cellFacei
1044  (
1045  IOobject
1046  (
1047  "cellFace",
1048  time().name(),
1049  *this,
1052  )
1053  );
1054 
1055  if (!cellFacei.empty())
1056  {
1057  checkFieldSize(cellFacei);
1058  facei_.resize(cellFacei.size());
1059  forAll(facei_, i)
1060  {
1061  facei_[i] = mesh_.cells()[celli_[i]][cellFacei[i]];
1062  }
1063  }
1064  }
1065 
1066  checkFieldSize(facei_);
1067  checkFieldSize(faceTrii_);
1068 
1069  // Request the tet base points so that they are built on all processors.
1070  // Constructing tet base points requires communication, so we can't leave
1071  // it until the tracking requests them as those calls are not synchronised.
1072  // Some processors might not be doing any tracking at all.
1073  mesh_.tetBasePtIs();
1074 
1075  // Mark the need to store the old-time cell-centres if the mesh is moving
1076  mesh_.oldCellCentres();
1077 }
1078 
1079 
1081 (
1083  const LagrangianState state
1084 )
1085 :
1086  mesh_(mesh)
1087 {
1088  mesh_.statesPtr_.set
1089  (
1091  );
1092 
1093  mesh_.offsetsPtr_.set
1094  (
1095  new labelList(mesh_.nGroups() + 1, 0)
1096  );
1097 
1098  for
1099  (
1100  label groupi = mesh_.stateToGroupi(state) + 1;
1101  groupi < mesh_.nGroups() + 1;
1102  ++ groupi
1103  )
1104  {
1105  mesh_.offsetsPtr_()[groupi] = mesh_.size();
1106  }
1107 
1108  constructNonConformal();
1109 
1110  mesh_.subMeshIndex_ = 0;
1111 
1112  Info<< indent;
1113  mesh_.printGroups(true);
1114  Info<< endl << indent;
1115  mesh_.printGroups(false);
1116  Info<< endl;
1117 
1118  forAll(mesh_.boundary(), patchi)
1119  {
1120  mesh_.boundary()[patchi].partition();
1121  }
1122 
1123  constructBehind();
1124 }
1125 
1126 
1128 (
1131 )
1132 :
1133  mesh_(mesh)
1134 {
1135  mesh_.statesPtr_.set
1136  (
1138  );
1139 
1140  mesh_.offsetsPtr_.set
1141  (
1142  new labelList(mesh_.nGroups() + 1, 0)
1143  );
1144 
1145  constructNonConformal();
1146 
1147  mesh_.subMeshIndex_ = 0;
1148 
1149  Info<< indent;
1150  mesh_.printGroups(true);
1151  Info<< endl;
1152  mesh_.partition();
1153 
1154  constructBehind();
1155 }
1156 
1157 
1158 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1159 
1161 {}
1162 
1163 
1165 {
1166  mesh_.statesPtr_.clear();
1167  mesh_.offsetsPtr_.clear();
1168 
1169  mesh_.origPatchNccPatchisPtr_.clear();
1170  mesh_.origPatchNccPatchesPtr_.clear();
1171  mesh_.nccPatchProcNccPatchisPtr_.clear();
1172  mesh_.receivePatchFacePtr_.clear();
1173  mesh_.receivePositionPtr_.clear();
1174 
1175  mesh_.fractionBehindPtr_.clear();
1176  mesh_.nTracksBehindPtr_.clear();
1177 
1178  mesh_.subMeshIndex_ = 0;
1179 
1180  forAll(mesh_.boundary(), patchi)
1181  {
1182  mesh_.boundary()[patchi].partition();
1183  }
1184 }
1185 
1186 
1188 (
1190 )
1191 :
1192  linear_(linear)
1193 {}
1194 
1195 
1197 (
1199  const LagrangianSubVectorField& quadratic
1200 )
1201 :
1202  linear_(linear),
1203  quadratic_(quadratic)
1204 {}
1205 
1206 
1207 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
1208 
1210 {
1211  if (!schemesPtr_.valid())
1212  {
1213  schemesPtr_ = new LagrangianSchemes(*this);
1214  }
1215 
1216  return schemesPtr_();
1217 }
1218 
1219 
1221 {
1222  if (!solutionPtr_.valid())
1223  {
1224  solutionPtr_ = new LagrangianSolution(*this);
1225  }
1226 
1227  return solutionPtr_();
1228 }
1229 
1230 
1232 {
1233  // Determine the number of global groups. Assumes processor patches come
1234  // after all global patches.
1235  label nGlobalGroups = static_cast<label>(LagrangianGroup::onPatchZero);
1236  forAll(boundary(), patchi)
1237  {
1238  const polyPatch& pp = boundary()[patchi].patch();
1239  if (isA<processorPolyPatch>(pp)) break;
1240  nGlobalGroups ++;
1241  }
1242 
1243  // Extend by one for the to-be-removed group
1244  nGlobalGroups ++;
1245 
1246  static const label completei =
1247  static_cast<label>(LagrangianGroup::complete);
1248  static const label inInternalMeshi =
1249  static_cast<label>(LagrangianGroup::inInternalMesh);
1250  static const label onPatchZeroi =
1251  static_cast<label>(LagrangianGroup::onPatchZero);
1252 
1253  // Create a list of sizes of the global groups
1254  labelList sizes(nGlobalGroups, 0);
1255  sizes[completei] = sub(LagrangianGroup::complete).size();
1256  sizes[inInternalMeshi] = sub(LagrangianGroup::inInternalMesh).size();
1257  forAll(boundary(), patchi)
1258  {
1259  const polyPatch& pp = boundary()[patchi].patch();
1260  if (isA<processorPolyPatch>(pp)) break;
1261  sizes[onPatchZeroi + patchi] = boundary()[patchi].mesh().size();
1262  }
1263  // (to-be-removed group last)
1264  sizes.last() =
1265  size() - boundary()[nGlobalGroups - 2 - onPatchZeroi].mesh().start();
1266 
1267  // Sum over all processes
1270 
1271  // Construct the result. This includes processor patch groups, on which a
1272  // "size" of -1 is set.
1273  labelList result(nGroups(), -1);
1274  for (label i = 0; i < nGlobalGroups - 1; ++ i)
1275  {
1276  result[i] = sizes[i];
1277  }
1278  // (to-be-removed group last)
1279  result.last() = sizes.last();
1280 
1281  return result;
1282 }
1283 
1284 
1287 {
1290  (
1291  positionName,
1292  *this,
1293  dimLength
1294  );
1295  LagrangianVectorInternalField& result = tresult.ref();
1296 
1297  forAll(coordinates_, i)
1298  {
1299  result[i] = position(i);
1300  }
1301 
1302  return tresult;
1303 }
1304 
1305 
1307 {
1308  return
1310  (
1311  mesh_,
1312  coordinates_[i], celli_[i], facei_[i], faceTrii_[i], 1
1313  );
1314 }
1315 
1316 
1318 (
1319  const point& position,
1321  label& celli,
1322  label& facei,
1323  label& faceTrii,
1324  const scalar fraction
1325 ) const
1326 {
1327  // Look for a containing cell and set the process if found
1328  remote procCelli;
1329  procCelli.elementi = mesh().cellTree().findInside(position);
1330  procCelli.proci = procCelli.elementi >= 0 ? Pstream::myProcNo() : -1;
1331 
1332  // Pick a unique processor
1333  reduce(procCelli, remote::firstProcOp());
1334 
1335  // Find the tetrahedron and the local coordinates
1337  if (procCelli.proci == Pstream::myProcNo())
1338  {
1339  result =
1341  (
1342  mesh_, position,
1343  coordinates, celli, facei, faceTrii, fraction
1344  )
1347  }
1348 
1349  // Communicate the location flag and return
1351 
1352  return result;
1353 }
1354 
1355 
1357 (
1358  const List<point>& position,
1360  labelList& celli,
1361  labelList& facei,
1363  const scalarList& fraction
1364 ) const
1365 {
1366  // Look for containing cells and set the process if found
1367  List<remote> procCelli(position.size());
1368  forAll(position, i)
1369  {
1370  procCelli[i].elementi = mesh().cellTree().findInside(position[0]);
1371  procCelli[i].proci =
1372  procCelli[i].elementi >= 0 ? Pstream::myProcNo() : -1;
1373  }
1374 
1375  // Pick unique processors
1376  reduce(procCelli, ListOp<remote::firstProcOp>());
1377 
1378  // Find the tetrahedra and the local coordinates
1380  forAll(position, i)
1381  {
1382  if (procCelli[i].proci != Pstream::myProcNo()) continue;
1383 
1384  result =
1386  (
1387  mesh_, position[i],
1388  coordinates[i], celli[i], facei[i], faceTrii[i], fraction[i]
1389  )
1392  }
1393 
1394  // Communicate the location flags and return
1396 
1397  return result;
1398 }
1399 
1400 
1402 {
1403  clearPosition();
1404 
1405  // Get the offsets
1406  checkPtr(offsetsPtr_, "Offsets");
1407  labelList& offsets = offsetsPtr_();
1408 
1409  // Partition the state offsets and states and create a permutation
1410  labelList permutation;
1411  switch (partitioningAlgorithm_)
1412  {
1414  permutation = partitionBin(offsets, states());
1415  break;
1416 
1418  permutation = partitionQuick(offsets, states());
1419  break;
1420  }
1421 
1422  // Print the updated states
1423  Info<< indent;
1424  printGroups(false);
1425  Info<< endl;
1426 
1427  // Apply the permutation to the states and positions
1428  permuteList(permutation, states());
1429  permuteList(permutation, coordinates_);
1430  permuteList(permutation, celli_);
1431  permuteList(permutation, facei_);
1432  permuteList(permutation, faceTrii_);
1433 
1434  // Apply the permutation to the non-conformal receive information (if any)
1435  if (receivePatchFacePtr_.valid())
1436  {
1437  permuteList(permutation, receivePatchFacePtr_());
1438  }
1439  if (receivePositionPtr_.valid())
1440  {
1441  permuteList(permutation, receivePositionPtr_());
1442  }
1443 
1444  // Check that the states are ordered
1445  #ifdef FULLDEBUG
1446  for (label groupi = 0; groupi < nGroups(); ++ groupi)
1447  {
1448  for (label i = offsets[groupi]; i < offsets[groupi + 1]; ++ i)
1449  {
1450  if (stateToGroupi(states()[i]) != groupi)
1451  {
1453  << "Partitioning failed"
1454  << exit(FatalError);
1455  }
1456  }
1457  }
1458  #endif
1459 
1460  // Reclaim space by removing the toBeRemoved group
1461  offsets[nGroups()] = offsets[nGroups() - 1];
1462 
1463  // Resize the states and positions
1464  resizeContainer(states());
1465  resizeContainer(coordinates_);
1466  resizeContainer(celli_);
1467  resizeContainer(facei_);
1468  resizeContainer(faceTrii_);
1469 
1470  // Resize the non-conformal receive information (if any)
1471  if (receivePatchFacePtr_.valid())
1472  {
1473  resizeContainer(receivePatchFacePtr_());
1474  }
1475  if (receivePositionPtr_.valid())
1476  {
1477  resizeContainer(receivePositionPtr_());
1478  }
1479 
1480  // Permute and resize the fields
1481  permuteAndResizeFields(permutation);
1482 
1483  // Update the patches
1484  forAll(boundary(), patchi)
1485  {
1486  boundary()[patchi].partition();
1487  }
1488 
1489  // Update the sub-all mesh
1490  subAll_.size_ = size();
1491 }
1492 
1493 
1494 template<class Displacement>
1496 (
1497  const List<LagrangianState>& endState,
1498  const Displacement& displacement,
1499  const LagrangianSubScalarField& deltaFraction,
1500  LagrangianSubScalarSubField& fraction
1501 )
1502 {
1503  clearPosition();
1504 
1505  // The fraction is about to change. Ensure the previous values are stored
1506  // to facilitate subsequent calculations.
1507  fraction.oldTime();
1508 
1509  // Track each element in the sub-mesh in turn
1510  forAll(fraction, subi)
1511  {
1512  const label i = subi + fraction.mesh().start();
1513 
1514  // Track to completion or the next face
1515  Tuple2<bool, scalar> onFaceAndF =
1517  (
1518  mesh_, displacement(subi), deltaFraction[subi],
1519  coordinates_[i], celli_[i], facei_[i], faceTrii_[i],
1520  fraction[subi],
1521  fractionBehindPtr_()[i], nTracksBehindPtr_()[i],
1522  debug
1523  ? static_cast<const string&>(name() + " #" + Foam::name(i))
1524  : NullObjectRef<string>()
1525  );
1526 
1527  // Update the state
1528  if (!onFaceAndF.first())
1529  {
1530  states()[i] = endState[subi];
1531  }
1532  else if (mesh_.isInternalFace(facei_[i]))
1533  {
1535  }
1536  else // if (<on a boundary face>)
1537  {
1538  // Determine the index of the patch that was tracked to
1539  label patchi =
1540  mesh_.boundaryMesh().patchIndices()
1541  [
1542  facei_[i] - mesh_.nInternalFaces()
1543  ];
1544 
1545  // If this patch has non-conformal cyclics associated with it, then
1546  // search through them and see if any was hit. If we find one that
1547  // does, override the patch index variable.
1548  if
1549  (
1550  origPatchNccPatchisPtr_.valid()
1551  && origPatchNccPatchisPtr_()[patchi].size()
1552  )
1553  {
1554  // Get the current position
1555  const point sendPosition =
1557  (
1558  mesh_,
1559  coordinates_[i], celli_[i], facei_[i], faceTrii_[i],
1560  fraction[subi]
1561  );
1562 
1563  // Get the displacement of the location that was hit
1564  const vector sendDisplacement =
1566  (
1567  mesh_,
1568  coordinates_[i], celli_[i], facei_[i], faceTrii_[i],
1569  fraction[subi]
1570  ).second();
1571 
1572  // Use ray searching on each non-conformal cyclic in turn
1573  forAll(origPatchNccPatchisPtr_()[patchi], patchNccPatchi)
1574  {
1575  const label nccPatchi =
1576  origPatchNccPatchisPtr_()[patchi][patchNccPatchi];
1577  const nonConformalCyclicPolyPatch& nccPp =
1578  origPatchNccPatchesPtr_()[patchi][patchNccPatchi];
1579 
1580  point receivePosition;
1581  const remote receiveProcAndFace =
1582  nccPp.ray
1583  (
1584  fraction[subi],
1585  nccPp.origPatch().whichFace(facei_[i]),
1586  sendPosition,
1587  displacement(subi, onFaceAndF.second())
1588  - fraction[subi]*sendDisplacement,
1589  receivePosition
1590  );
1591 
1592  const label receiveProci = receiveProcAndFace.proci;
1593 
1594  if (receiveProci == -1) continue;
1595 
1596  const label receiveFacei = receiveProcAndFace.elementi;
1597 
1598  receivePatchFacePtr_()[i] = receiveFacei;
1599  receivePositionPtr_()[i] = receivePosition;
1600 
1601  patchi =
1602  nccPatchProcNccPatchisPtr_()[nccPatchi][receiveProci];
1603 
1604  break;
1605  }
1606  }
1607 
1608  // Set the state to that of the identified patch
1609  states()[i] =
1610  static_cast<LagrangianState>
1611  (
1612  static_cast<label>(LagrangianState::onPatchZero)
1613  + patchi
1614  );
1615  }
1616  }
1617 }
1618 
1619 
1620 template
1621 void Foam::LagrangianMesh::track<Foam::LagrangianMesh::linearDisplacement>
1622 (
1623  const List<LagrangianState>& endState,
1624  const linearDisplacement& displacement,
1625  const LagrangianSubScalarField& deltaFraction,
1626  LagrangianSubScalarSubField& fraction
1627 );
1628 
1629 
1630 template
1631 void Foam::LagrangianMesh::track<Foam::LagrangianMesh::parabolicDisplacement>
1632 (
1633  const List<LagrangianState>& endState,
1634  const parabolicDisplacement& displacement,
1635  const LagrangianSubScalarField& deltaFraction,
1636  LagrangianSubScalarSubField& fraction
1637 );
1638 
1639 
1641 (
1642  const LagrangianScalarInternalDynamicField& fraction
1643 )
1644 {
1645  clearPosition();
1646 
1647  // Internal face-crossings
1648  const LagrangianSubMesh incompleteMesh
1649  (
1651  );
1652 
1653  forAll(incompleteMesh, subi)
1654  {
1655  const label i = subi + incompleteMesh.start();
1656 
1657  if (states()[i] != LagrangianState::onInternalFace) continue;
1658 
1659  // Cross the face
1661  (
1662  mesh_,
1663  coordinates_[i], celli_[i], facei_[i], faceTrii_[i]
1664  );
1665 
1666  // Update the state
1668  }
1669 
1670  // Patch-face crossings and boundary condition evaluations
1671  if
1672  (
1675  )
1676  {
1678 
1679  // Initialise patch crossings
1680  forAll(boundary(), patchi)
1681  {
1682  boundary()[patchi].initEvaluate(pBufs, *this, fraction);
1683  }
1684 
1685  // Initialise patch field evaluations
1686  #define INIT_EVAL_TYPE_PATCH_FIELDS(Type, GeoField) \
1687  { \
1688  HashTable<GeoField<Type>*> fields \
1689  ( \
1690  lookupCurrentFields<GeoField<Type>>() \
1691  ); \
1692  \
1693  forAllIter(typename HashTable<GeoField<Type>*>, fields, iter) \
1694  { \
1695  forAll(boundary(), patchi) \
1696  { \
1697  iter()->boundaryFieldRef()[patchi].initEvaluate \
1698  ( \
1699  pBufs, \
1700  fraction \
1701  ); \
1702  } \
1703  } \
1704  }
1709  #undef INIT_EVAL_TYPE_PATCH_FIELDS
1710 
1711  // Block for any outstanding requests
1712  pBufs.finishedSends();
1713 
1714  // Patch crossings
1715  forAll(boundary(), patchi)
1716  {
1717  boundary()[patchi].evaluate(pBufs, *this, fraction);
1718  }
1719 
1720  // Patch field evaluations
1721  #define EVAL_TYPE_PATCH_FIELDS(Type, GeoField) \
1722  { \
1723  HashTable<GeoField<Type>*> fields \
1724  ( \
1725  lookupCurrentFields<GeoField<Type>>() \
1726  ); \
1727  \
1728  forAllIter(typename HashTable<GeoField<Type>*>, fields, iter) \
1729  { \
1730  forAll(boundary(), patchi) \
1731  { \
1732  iter()->boundaryFieldRef()[patchi].evaluate \
1733  ( \
1734  pBufs, \
1735  fraction \
1736  ); \
1737  } \
1738  } \
1739  }
1744  #undef EVAL_TYPE_PATCH_FIELDS
1745  }
1746  else
1747  {
1749  << "Unsupported communications type "
1751  << exit(FatalError);
1752  }
1753 }
1754 
1755 
1756 void Foam::LagrangianMesh::reset(const bool initial, const bool final)
1757 {
1758  clearPosition();
1759 
1760  if (initial && final) return;
1761 
1762  if (initial)
1763  {
1764  coordinates_.storeOldTimes();
1765  coordinates_.oldTime();
1766  celli_.storeOldTimes();
1767  celli_.oldTime();
1768  facei_.storeOldTimes();
1769  facei_.oldTime();
1770  faceTrii_.storeOldTimes();
1771  faceTrii_.oldTime();
1772 
1773  #define OLD_TIME_TYPE_FIELDS(Type, GeoField) \
1774  { \
1775  HashTable<GeoField<Type>*> fields \
1776  ( \
1777  lookupCurrentFields<GeoField<Type>>() \
1778  ); \
1779  \
1780  forAllIter(typename HashTable<GeoField<Type>*>, fields, iter) \
1781  { \
1782  iter()->storeOldTimes(); \
1783  iter()->oldTime(); \
1784  } \
1785  }
1788  #undef OLD_TIME_TYPE_FIELDS
1789  }
1790 
1791  if (!initial)
1792  {
1793  static_cast<DynamicField<barycentric>&>(coordinates_) =
1794  coordinates_.oldTime();
1795  static_cast<DynamicField<label>&>(celli_) = celli_.oldTime();
1796  static_cast<DynamicField<label>&>(facei_) = facei_.oldTime();
1797  static_cast<DynamicField<label>&>(faceTrii_) = faceTrii_.oldTime();
1798 
1799  #define RESET_OLD_TIME_TYPE_FIELDS(Type, GeoField) \
1800  { \
1801  HashTable<GeoField<Type>*> fields \
1802  ( \
1803  lookupCurrentFields<GeoField<Type>>() \
1804  ); \
1805  \
1806  forAllIter(typename HashTable<GeoField<Type>*>, fields, iter) \
1807  { \
1808  iter()->reset(iter()->oldTime()); \
1809  } \
1810  }
1813  #undef RESET_OLD_TIME_TYPE_FIELDS
1814  }
1815 
1816  if (final)
1817  {
1818  coordinates_.clearOldTimes();
1819  celli_.clearOldTimes();
1820  facei_.clearOldTimes();
1821  faceTrii_.clearOldTimes();
1822 
1823  #define CLEAR_OLD_TIME_TYPE_FIELDS(Type, GeoField) \
1824  { \
1825  HashTable<GeoField<Type>*> fields \
1826  ( \
1827  lookupCurrentFields<GeoField<Type>>() \
1828  ); \
1829  \
1830  forAllIter(typename HashTable<GeoField<Type>*>, fields, iter) \
1831  { \
1832  iter()->clearOldTimes(); \
1833  } \
1834  }
1837  #undef CLEAR_OLD_TIME_TYPE_FIELDS
1838  }
1839 
1841  #define INSERT_FIELD_NAMES(Type, GeoField) \
1842  fieldNames.insert(lookupCurrentFields<GeoField<Type>>(true).toc());
1849  #undef INSERT_FIELD_NAMES
1850 
1851  if (!fieldNames.empty())
1852  {
1854  << "Non-dynamic and/or internal fields " << fieldNames.sortedToc()
1855  << " are registered during a reset. Only fields relating to "
1856  << "fundamental state should be present at this time, and these "
1857  << "should all be dynamic non-internal fields" << exit(FatalError);
1858  }
1859 
1860  subAll_.size_ = size();
1861 }
1862 
1863 
1865 {
1866  clearPosition();
1867 
1868  if (statesPtr_.valid())
1869  {
1870  states().clear();
1871  }
1872 
1873  coordinates_.clear();
1874  celli_.clear();
1875  facei_.clear();
1876  faceTrii_.clear();
1877 
1878  subAll_.size_ = 0;
1879 }
1880 
1881 
1883 (
1884  const label nGroups,
1885  const UList<labelPair>& elementsGroups
1886 )
1887 {
1888  if (elementsGroups.empty()) return labelList(nGroups + 1, size());
1889 
1890  // !!! Not implemented as yet. This is needed for models which
1891  // instantaneously modify or remove elements. The only instantaneous models
1892  // implemented thus far are injection models which only create elements.
1893  // This will need implementing for breakup models and similar.
1895  return labelList();
1896 }
1897 
1898 
1900 {
1901  List<labelPair> elementsGroups(elements.size());
1902  forAll(elements, i)
1903  {
1904  elementsGroups[i] = labelPair(elements[i], 0);
1905  }
1906 
1907  const labelList offsets = partition(1, elementsGroups);
1908 
1909  remove(offsets[1] - offsets[0]);
1910 }
1911 
1912 
1914 {
1915  if (nElements == 0) return;
1916 
1917  // !!! Not implemented as yet. This is needed for models which
1918  // instantaneously remove elements. The only instantaneous models
1919  // implemented thus far are injection models which only create elements.
1920  // This will need implementing for breakup models and similar.
1922 }
1923 
1924 
1926 {
1927  if (foundObject<LagrangianVectorInternalField>(positionName))
1928  {
1929  lookupObjectRef<LagrangianVectorInternalField>(positionName).checkOut();
1930  }
1931 }
1932 
1933 
1935 {
1936  if (!this->foundObject<LagrangianVectorInternalField>(positionName))
1937  {
1938  this->position().ptr()->store();
1939  }
1940 }
1941 
1942 
1944 {
1945  if (map.reverseCellMap().empty()) return;
1946 
1948 
1949  meshObjects::topoChange<LagrangianMesh>(*this, map);
1950 }
1951 
1952 
1954 {
1956 
1957  meshObjects::mapMesh<LagrangianMesh>(*this, map);
1958 }
1959 
1960 
1962 {
1964 
1965  meshObjects::distribute<LagrangianMesh>(*this, map);
1966 }
1967 
1968 
1970 (
1974  const bool write
1975 ) const
1976 {
1977  return objectRegistry::writeObject(fmt, ver, cmp, write);
1978 }
1979 
1980 
1981 bool Foam::LagrangianMesh::write(const bool write) const
1982 {
1983  return objectRegistry::write(write);
1984 }
1985 
1986 
1987 // ************************************************************************* //
Istream and Ostream manipulators taking arguments.
#define RESET_OLD_TIME_TYPE_FIELDS(Type, GeoField)
#define INIT_EVAL_TYPE_PATCH_FIELDS(Type, GeoField)
#define INJECT_TYPE_FIELDS(Type, GeoField)
#define CLEAR_OLD_TIME_TYPE_FIELDS(Type, GeoField)
#define BIRTH_TYPE_FIELDS(Type, GeoField)
#define OLD_TIME_TYPE_FIELDS(Type, GeoField)
#define PERMUTE_TYPE_FIELDS(Type, GeoField)
#define INJECT_STATE_FIELD(GeoField)
#define INSERT_INTERNAL_FIELD_NAMES(Type, GeoField)
#define INSERT_FIELD_NAMES(Type, GeoField)
#define EVAL_TYPE_PATCH_FIELDS(Type, GeoField)
Various functions to operate on Lists.
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const
Return mesh.
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Dynamically sized Field.
Definition: DynamicField.H:72
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
Generic GeometricField class.
A HashTable with keys but without contents.
Definition: HashSet.H:62
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
readOption
Enumeration defining the read options.
Definition: IOobject.H:117
writeOption writeOpt() const
Definition: IOobject.H:367
const word & name() const
Return name.
Definition: IOobject.H:307
writeOption
Enumeration defining the write options.
Definition: IOobject.H:126
Version number type.
Definition: IOstream.H:97
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:87
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:194
const LagrangianMesh & mesh() const
Return the mesh reference.
changer(LagrangianMesh &mesh, const LagrangianState state)
Construct for a Lagrangian mesh with a given state.
linearDisplacement(const LagrangianSubVectorField &linear)
Construct from a reference to the displacement.
parabolicDisplacement(const LagrangianSubVectorField &linear, const LagrangianSubVectorField &quadratic)
Construct from references to the displacements.
Class containing Lagrangian geometry and topology.
LagrangianState state(const label i) const
Return the state for an element of the mesh, or a none state.
const Time & time() const
Return time.
LagrangianMesh(const polyMesh &mesh, const word &name, const IOobject::readOption readOption=IOobject::READ_IF_PRESENT, const IOobject::writeOption writeOption=IOobject::AUTO_WRITE)
Construct from a mesh and a name.
const labelIODynamicField & faceTrii() const
Access the face-tet indices.
const polyMesh & mesh() const
Access the mesh.
label size() const
Get the number of elements.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
location
Enumeration for the locations of searched positions.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
static const NamedEnum< permutationAlgorithm, 2 > permutationAlgorithmNames_
Permutation algorithm names.
~LagrangianMesh()
Destructor.
static permutationAlgorithm permutationAlgorithm_
Permutation algorithm.
static const word stateName
Name of the state field.
virtual void storePosition()
Store the positions for use during mapping.
partitioningAlgorithm
Enumeration of the partitioning algorithm.
void track(const List< LagrangianState > &endState, const Displacement &displacement, const LagrangianSubScalarField &deltaFraction, LagrangianSubScalarSubField &fraction)
Track the positions along the given displacements.
const labelIODynamicField & celli() const
Access the cell indices.
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write=true) const
Write using given format, version and compression.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
location locate(const point &position, barycentric &coordinates, label &celli, label &facei, label &faceTrii, const scalar fraction) const
Convert a position into a set of coordinates and a.
static const word coordinatesName
Name of the coordinates field.
virtual void clearPosition()
Clear the positions uses during mapping.
const LagrangianBoundaryMesh & boundary() const
Return reference to boundary mesh.
label nGroups() const
Return the number of groups.
virtual bool write(const bool write=true) const
Write settings from the database.
tmp< LagrangianVectorInternalField > position() const
Return the global positions.
static const word fractionName
Name of the tracked fraction field.
static const word positionName
Name of the position field.
void reset(const bool initial, const bool final)
Reset the mesh to the old-time conditions.
labelList subMeshGlobalSizes() const
Return the global sizes of all the sub-meshes. A value of -1.
const LagrangianSolution & solution() const
Access the solution controls.
permutationAlgorithm
Enumeration of the permutation algorithm.
void clear()
Clear all geometry out of the Lagrangian mesh.
const LagrangianSchemes & schemes() const
Access the schemes.
void crossFaces(const LagrangianScalarInternalDynamicField &fraction)
Cross the faces.
label stateToGroupi(const LagrangianState state) const
Convert a state to a group index.
static partitioningAlgorithm partitioningAlgorithm_
Partitioning algorithm.
void remove(const UList< label > &elements)
Remove specified elements from the mesh. Shuffles everything.
static const NamedEnum< partitioningAlgorithm, 2 > partitioningAlgorithmNames_
Partitioning algorithm names.
void partition()
Partition the mesh such that the groups are contiguous in memory.
LagrangianSubMesh sub(const LagrangianGroup group) const
Return a sub-mesh for the given group.
static const word prefix
Instance prefix.
const List< LagrangianState > & states() const
Return the states.
const barycentricIODynamicField & coordinates() const
Access the coordinates.
const labelIODynamicField & facei() const
Access the cell-face indices.
Selector class for Lagrangian schemes.
Selector class for Lagrangian schemes.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
label size() const
Return size.
label start() const
Return start.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static const List< word > & null()
Return a null List.
Definition: ListI.H:118
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
const Field0Type & oldTime() const
Return the old-time field.
Definition: OldTimeField.C:322
void storeOldTimes() const
Store the old-time fields.
Definition: OldTimeField.C:257
void clearOldTimes()
Clear old-time fields.
Definition: OldTimeField.C:276
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
const Type2 & second() const
Return second.
Definition: Tuple2.H:131
const Type1 & first() const
Return first.
Definition: Tuple2.H:119
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:74
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
T & last()
Return the last element of the list.
Definition: UListI.H:128
static const NamedEnum< commsTypes, 3 > commsTypeNames
Definition: UPstream.H:71
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static const Form nan
Definition: VectorSpace.H:124
static direction size()
Return the number of elements in the VectorSpace = Ncmpts.
Definition: VectorSpaceI.H:83
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
Centred interpolation interpolation scheme class.
Definition: linear.H:53
Non-conformal cyclic poly patch. As nonConformalCoupledPolyPatch, but the neighbouring patch is local...
remote ray(const scalar fraction, const label origFacei, const vector &p, const vector &n, point &nbrP) const
Compute a ray intersection across the coupling.
const polyPatch & origPatch() const
Original patch.
Registry of regIOobjects.
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write) const
Write the objects.
void clear()
Remove all regIOobject owned by the registry.
const labelList & patchIndices() const
Boundary face patch indices.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:1046
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:1080
virtual const pointField & oldCellCentres() const
Return old cell centres for mesh motion.
Definition: polyMesh.C:1387
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:360
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reverseCellMap() const
Reverse cell map.
label nInternalFaces() const
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
const cellList & cells() const
virtual bool write(const bool write=true) const
Write using setting from DB.
Struct for keeping processor, element (cell, face, point) index.
Definition: remote.H:57
label elementi
Element index.
Definition: remote.H:70
label proci
Processor index.
Definition: remote.H:67
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
static List< word > fieldNames
Definition: globalFoam.H:46
Enum namedEnumOptimisationSwitch(const char *name, const NamedEnum< Enum, nEnum > &enumNames, const Enum defaultValue)
Lookup optimisation switch or add default value.
Definition: debug.H:94
barycentric coordinates(const polyMesh &mesh, const point &position, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the coordinates given the position and tet topology.
Definition: tracking.C:1259
Pair< vector > faceNormalAndDisplacement(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the normal of the corresponding point on the associated face and.
Definition: tracking.C:1302
void crossInternalFace(const polyMesh &mesh, barycentric &coordinates, label &celli, label &facei, label &faceTrii)
Cross an internal face.
Definition: tracking.C:1681
point position(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the position given the coordinates and tet topology.
Definition: trackingI.H:224
Tuple2< bool, scalar > toFace(const polyMesh &mesh, const Displacement &displacement, const scalar fraction, barycentric &coordinates, label &celli, label &facei, label &faceTrii, scalar &stepFraction, scalar &stepFractionBehind, label &nTracksBehind, const string &debugPrefix=NullObjectRef< string >())
Track along the displacement for a given fraction of the overall.
bool locate(const polyMesh &mesh, const point &position, barycentric &coordinates, label &celli, label &facei, label &faceTrii, const scalar stepFraction, const string &debugPrefix=NullObjectRef< string >())
Initialise the location at the given position. Returns whether or not a.
Definition: tracking.C:1593
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
List< label > labelList
A List of labels.
Definition: labelList.H:56
Field< barycentric > barycentricField
Barycentric field.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
labelList invert(const label len, const labelUList &)
Invert one-to-one map. Unmapped elements will be -1.
Definition: ListOps.C:37
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
labelList findIndices(const ListType &, typename ListType::const_reference, const label start=0)
Find all occurrences of given element. Linear search.
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
LagrangianState
Lagrangian state enumeration.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
DimensionedField< Type, LagrangianMesh > LagrangianInternalField
defineTypeNameAndDebug(combustionModel, 0)
GeometricField< Type, LagrangianMesh, LagrangianPrimitiveDynamicField > LagrangianDynamicField
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
LagrangianSubSubField< scalar > LagrangianSubScalarSubField
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:210
LagrangianGroup
Lagrangian group enumeration.
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
void Swap(T &a, T &b)
Definition: Swap.H:43
LagrangianSubField< scalar > LagrangianSubScalarField
DimensionedField< Type, LagrangianMesh, LagrangianPrimitiveDynamicField > LagrangianInternalDynamicField
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
GeometricField< Type, LagrangianMesh > LagrangianField
volScalarField & p
Operator to apply a binary operation to a pair of lists.
Definition: ListOps.H:288
Operator to take the first valid process.
Definition: remote.H:93