LagrangianMesh.H
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 Class
25  Foam::LagrangianMesh
26 
27 Description
28  Class containing Lagrangian geometry and topology
29 
30 SourceFiles
31  LagrangianMesh.C
32  LagrangianMeshTemplates.C
33  LagrangianMeshI.H
34  LagrangianMeshLocation.H
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef LagrangianMesh_H
39 #define LagrangianMesh_H
40 
41 #include "GeoMesh.H"
42 #include "barycentricField.H"
44 #include "LagrangianBoundaryMesh.H"
45 #include "LagrangianFieldsFwd.H"
46 #include "LagrangianSchemes.H"
47 #include "LagrangianSolution.H"
48 #include "LagrangianState.H"
49 #include "LagrangianSubMesh.H"
50 #include "labelIODynamicField.H"
51 #include "polyMesh.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 namespace Foam
56 {
57 
58 class LagrangianInjection;
59 class nonConformalCyclicPolyPatch;
60 
61 /*---------------------------------------------------------------------------*\
62  Class LagrangianMesh Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 class LagrangianMesh
66 :
67  public objectRegistry,
68  public GeoMesh<polyMesh>
69 {
70 public:
71 
72  // Public Enumerations
73 
74  //- Enumeration of the permutation algorithm
75  enum class permutationAlgorithm
76  {
77  copy,
78  inPlace
79  };
80 
81  //- Permutation algorithm names
84 
85  //- Enumeration of the partitioning algorithm
86  enum class partitioningAlgorithm
87  {
88  bin,
89  quick
90  };
91 
92  //- Partitioning algorithm names
95 
96  //- Enumeration for the locations of searched positions
97  enum class location
98  {
99  inCell,
100  onBoundary,
102  };
103 
104 
105 private:
106 
107  // Private Data
108 
109  //- Reference to the mesh
110  const polyMesh& mesh_;
111 
112  //- Local coordinates of the elements within their tetrahedra
113  barycentricIODynamicField coordinates_;
114 
115  //- Indices of the cells that contain the tetrahedra
116  labelIODynamicField celli_;
117 
118  //- Indices of the faces that form the bases of the tetrahedra
119  labelIODynamicField facei_;
120 
121  //- Indices of the face-triangles that form the bases of the tetrahedra
122  labelIODynamicField faceTrii_;
123 
124  //- Boundary mesh
125  LagrangianBoundaryMesh boundary_;
126 
127  //- A special sub-mesh which spans the entire mesh
128  LagrangianSubMesh subAll_;
129 
130  //- The elements' states
132 
133  //- Offsets defining ranges associated with the elements' groups
134  autoPtr<labelList> offsetsPtr_;
135 
136  //- For every patch, a list of associated non-conformal cyclic patch
137  // indices. Empty if there are no non-conformal cyclic patches.
138  autoPtr<labelListList> origPatchNccPatchisPtr_;
139 
140  //- For every patch, a list of associated non-conformal cyclic patch
141  // pointers. Empty if there are no non-conformal cyclic patches.
143  origPatchNccPatchesPtr_;
144 
145  //- For every non-conformal cyclic patch, the indices of the associated
146  // non-conformal processor-cyclic patches. Empty if there are no
147  // non-conformal cyclic patches.
148  autoPtr<labelListList> nccPatchProcNccPatchisPtr_;
149 
150  //- The receiving patch face for non-conformal cyclic transfers
151  autoPtr<DynamicList<label>> receivePatchFacePtr_;
152 
153  //- The receiving position for non-conformal cyclic transfers
154  autoPtr<DynamicList<point>> receivePositionPtr_;
155 
156  //- The step fraction less than the maximum reached so far. See
157  // tracking.H for details.
158  autoPtr<LagrangianDynamicField<scalar>> fractionBehindPtr_;
159 
160  //- The number of tracks carried out that ended in a step fraction less
161  // than the maximum reached so far. See tracking.H for details.
162  autoPtr<LagrangianDynamicField<label>> nTracksBehindPtr_;
163 
164  //- Sub-mesh index
165  mutable label subMeshIndex_;
166 
167  //- Schemes created on demand
168  mutable autoPtr<LagrangianSchemes> schemesPtr_;
169 
170  //- Solution controls created on demand
171  mutable autoPtr<LagrangianSolution> solutionPtr_;
172 
173 
174  // Private Member Functions
175 
176  // Checks
177 
178  //- Check that the field is the same size as the mesh
179  template<class Field>
180  inline void checkFieldSize(const Field& field) const;
181 
182  //- Check that a pointer is valid
183  template<class Type>
184  inline void checkPtr
185  (
186  const autoPtr<Type>& ptr,
187  const word& name
188  ) const;
189 
190  //- Print a row of the group table
191  void printGroups(const bool header) const;
192 
193 
194  // Partitioning and Permutation
195 
196  //- Partition the states using the bin-sort-based algorithm
197  labelList partitionBin
198  (
199  labelList& offsets,
201  ) const;
202 
203  //- Partition the states using a quick-sort-based algorithm
204  labelList partitionQuick
205  (
206  labelList& offsets,
208  ) const;
209 
210  //- Reorder and resize all registered fields using the given
211  // permutation
212  void permuteAndResizeFields(const labelList& permutation);
213 
214  //- Reorder a list with the given permutation
215  template<class Type>
216  static void permuteList
217  (
218  const labelList& permutation,
219  UList<Type>& list
220  );
221 
222  //- Reorder a list using the copy algorithm
223  template<class Type>
224  static void permuteListCopy
225  (
226  const labelList& permutation,
227  UList<Type>& list
228  );
229 
230  //- Reorder a list using the in-place algorithm
231  template<class Type>
232  static void permuteListInPlace
233  (
234  const labelList& permutation,
235  UList<Type>& list
236  );
237 
238  //- Resize a container to match the mesh
239  template<class Container>
240  void resizeContainer(Container& container) const;
241 
242 
243  // Addition
244 
245  //- Append specified elements in the mesh with the given geometry
246  // and topology. Does nothing to fields.
247  LagrangianSubMesh append
248  (
250  const labelField& celli,
251  const labelField& facei,
252  const labelField& faceTrii
253  );
254 
255  //- Birth specified elements in the mesh from the given list of
256  // parent indices. Does nothing to fields.
257  LagrangianSubMesh append(const labelList& parents);
258 
259  //- Expand and set values for explicitly specified fields following
260  // an injection or birth event
261  template<class Type, template<class> class GeoField>
262  void appendSpecifiedField
263  (
264  const LagrangianSubMesh& appendMesh,
265  GeoField<Type>& geoField,
266  const Field<Type>& field
267  ) const;
268 
269  //- Expand and set values for explicitly specified fields following
270  // an injection or birth event
271  template<class Type, template<class> class GeoField>
272  bool appendSpecifiedField
273  (
274  const LagrangianSubMesh& appendMesh,
275  const word& fieldName,
276  const Field<Type>& field
277  ) const;
278 
279  //- Expand and set values for explicitly specified fields following
280  // an injection or birth event
281  template<class Type, class ... FieldNamesAndFields>
282  wordHashSet appendSpecifiedFields
283  (
284  const LagrangianSubMesh& appendMesh,
285  const word& fieldName,
286  const Field<Type>& field,
287  const FieldNamesAndFields& ... fieldNamesAndFields
288  ) const;
289 
290  //- Termination clause for the above
291  wordHashSet appendSpecifiedFields
292  (
293  const LagrangianSubMesh& appendMesh
294  ) const;
295 
296  //- Expand and set values for unspecified fields following an
297  // injection
298  void injectUnspecifiedFields
299  (
300  const LagrangianSubMesh& injectionMesh,
301  const wordHashSet& specifiedFieldNames
302  );
303 
304  //- Expand and set values for unspecified fields following an
305  // injection
306  void injectUnspecifiedFields
307  (
308  const LagrangianInjection& injection,
309  const LagrangianSubMesh& injectionMesh,
310  const wordHashSet& specifiedFieldNames
311  );
312 
313  //- Expand and set values for unspecified fields following a birth
314  void birthUnspecifiedFields
315  (
316  const labelList& parents,
317  const LagrangianSubMesh& birthMesh,
318  const wordHashSet& specifiedFieldNames
319  );
320 
321 
322 public:
323 
324  // Friend classes
325 
326  //- The non-conformal cyclic patch needs to access the receiving
327  // information populated during tracking
329 
330  //- The processor patch communicates and adds elements to the receiving
331  // side, but field values are added later by the corresponding patch
332  // field. This is not permitted by the public interface in which
333  // elements and associated field values are required to be added
334  // simultaneously so that they are guaranteed to remain consistent.
335  // The requirements of the processor patch and patch fields are
336  // considered a special case, so access to the private append methods
337  // is permitted to these classes.
338  friend class processorLagrangianPatch;
339 
340  //- See above
341  template<class Type>
342  friend class processorLagrangianPatchField;
343 
344  //- See above
346 
347  //- See above
348  template<class Type>
350 
351 
352  // Public classes
353 
354  //- Class to define the scope of Lagrangian mesh state changes
355  class changer
356  {
357  // Private Data
358 
359  //- Reference to the Lagrangian mesh
360  LagrangianMesh& mesh_;
361 
362 
363  // Private Member Functions
364 
365  //- Construct the mesh's non-conformal data
366  void constructNonConformal() const;
367 
368  //- Construct the mesh's behind data
369  void constructBehind() const;
370 
371 
372  public:
373 
374  // Constructors
375 
376  //- Construct for a Lagrangian mesh with a given state
377  changer
378  (
380  const LagrangianState state
381  );
382 
383  //- Construct for a Lagrangian mesh with given states
384  changer
385  (
388  );
389 
390 
391  // Destructor
392  ~changer();
393  };
394 
395  //- Class to hold and index into the field references associated with a
396  // linear displacement
397  class linearDisplacement
398  {
399  // Private Member Data
400 
401  //- Reference to the linear displacement
402  const LagrangianSubVectorField& linear_;
403 
404  public:
405 
406  // Constructors
407 
408  //- Construct from a reference to the displacement
410 
411 
412  // Member Operators
413 
414  //- The displacement for a given index
415  inline const vector& operator()(const label i) const;
416 
417  //- The displacement for a given index and remaining fraction
418  inline const vector& operator()
419  (
420  const label i,
421  const scalar f
422  ) const;
423  };
424 
425  //- Class to hold and index into the field references associated with a
426  // parabolic displacement
428  {
429  // Private Member Data
430 
431  //- Reference to the linear displacement
432  const LagrangianSubVectorField& linear_;
433 
434  //- Reference to the quadratic displacement
435  const LagrangianSubVectorField& quadratic_;
436 
437 
438  public:
439 
440  // Constructors
441 
442  //- Construct from references to the displacements
444  (
446  const LagrangianSubVectorField& quadratic
447  );
448 
449 
450  // Member Operators
451 
452  //- The displacements for a given index
453  inline Pair<vector> operator()(const label i) const;
454 
455  //- The displacement for a given index and remaining fraction
456  inline vector operator()
457  (
458  const label i,
459  const scalar f
460  ) const;
461  };
462 
463  //- Class to hold element-group indices, and associate the group
464  // indices with a given enumeration. Binary compatible with labelPair
465  // so that a UList<elementGroup<Enumeration>> can be cast to a
466  // UList<labelPair>. This allows the labelPair-based partition
467  // function to be used with an arbitrary enumeration without
468  // templating the implementation.
469  template<class Enumeration>
470  class elementGroup
471  :
472  private labelPair
473  {
474  //- Reference class returned by the group accessors
475  class EnumerationRef
476  {
477  private:
478 
479  label& i_;
480 
481  public:
482 
483  EnumerationRef(label& i)
484  :
485  i_(i)
486  {}
487 
488  operator Enumeration() const
489  {
490  return static_cast<Enumeration>(i_);
491  }
492 
493  void operator=(const Enumeration e)
494  {
495  i_ = static_cast<label>(e);
496  }
497  };
498 
499  //- Modify the element index
500  label& element()
501  {
502  return first();
503  }
504 
505  //- Access the element index
506  label element() const
507  {
508  return first();
509  }
510 
511  //- Modify the group enumeration
512  EnumerationRef group()
513  {
514  return EnumerationRef(second());
515  }
516 
517  //- Access the group enumeration
518  Enumeration group() const
519  {
520  return static_cast<Enumeration>(second());
521  }
522  };
523 
524 
525  // Public Static Data
526 
527  //- Run-time type information
528  ClassName("LagrangianMesh");
529 
530  //- Instance prefix
531  static const word prefix;
532 
533  //- Name of the coordinates field
534  static const word coordinatesName;
535 
536  //- Name of the position field
537  static const word positionName;
538 
539  //- Name of the state field
540  static const word stateName;
541 
542  //- Name of the tracked fraction field
543  static const word fractionName;
544 
545  //- Permutation algorithm
547 
548  //- Partitioning algorithm
550 
551 
552  // Public Type Definitions
553 
554  //- Mesh type
555  typedef LagrangianMesh Mesh;
556 
557  //- Boundary mesh type
559 
560  //- Patch field type
561  template<class Type>
563 
564  //- Field source type
565  template<class Type>
567 
568 
569  // Constructors
570 
571  //- Construct from a mesh and a name
573  (
574  const polyMesh& mesh,
575  const word& name,
578  );
579 
580  //- Construct from a mesh and a name and a list of patch types
582  (
583  const polyMesh& mesh,
584  const word& name,
585  const wordList& wantedPatchTypes,
588  );
589 
590  //- Disallow default bitwise copy construction
591  LagrangianMesh(const LagrangianMesh&) = delete;
592 
593 
594  //- Destructor
595  ~LagrangianMesh();
596 
597 
598  // Member Functions
599 
600  // Access
601 
602  using objectRegistry::time;
604  using objectRegistry::operator!=;
605  using objectRegistry::operator==;
606 
607  //- Access the mesh
608  inline const polyMesh& mesh() const;
609 
610  //- Access the coordinates
611  inline const barycentricIODynamicField& coordinates() const;
612 
613  //- Access the cell indices
614  inline const labelIODynamicField& celli() const;
615 
616  //- Access the cell-face indices
617  inline const labelIODynamicField& facei() const;
618 
619  //- Access the face-tet indices
620  inline const labelIODynamicField& faceTrii() const;
621 
622  //- Return reference to boundary mesh
623  inline const LagrangianBoundaryMesh& boundary() const;
624 
625  //- Return whether or not the mesh is changing
626  inline bool changing() const;
627 
628  //- Return the states
629  inline const List<LagrangianState>& states() const;
630 
631  //- Access the states
632  inline List<LagrangianState>& states();
633 
634  //- Return the state for an element of the mesh, or a none state
635  // if not changing
636  inline LagrangianState state(const label i) const;
637 
638  //- Return the state for an element of a sub mesh, or a none state
639  // if not changing
640  inline LagrangianState state
641  (
642  const LagrangianSubMesh& subMesh,
643  const label subi
644  ) const;
645 
646  //- Get a sub-mesh index
647  inline label subMeshIndex() const;
648 
649  //- Access the schemes
650  const LagrangianSchemes& schemes() const;
651 
652  //- Access the solution controls
653  const LagrangianSolution& solution() const;
654 
655  //- Lookup all current-time fields of the given type
656  template<class GeoField>
658  (
659  const bool strict = false
660  ) const;
661 
662  //- Lookup all current-time fields of the given type
663  template<class GeoField>
665  (
666  const bool strict = false
667  );
668 
669 
670  // Check
671 
672  //- Get the number of elements
673  inline label size() const;
674 
675  //- Get the global number of elements
676  inline label globalSize() const;
677 
678  //- Return size
679  static inline label size(const LagrangianMesh& mesh);
680 
681  //- Return the number of states
682  inline label nStates() const;
683 
684  //- Return the number of groups
685  inline label nGroups() const;
686 
687  //- Convert a state to a group index
688  inline label stateToGroupi(const LagrangianState state) const;
689 
690  //- Return a sub-mesh for the given group
691  inline LagrangianSubMesh sub(const LagrangianGroup group) const;
692 
693  //- Return a sub-mesh for no elements
694  inline LagrangianSubMesh subNone() const;
695 
696  //- Return a sub-mesh for all elements
697  inline const LagrangianSubMesh& subAll() const;
698 
699  //- Return a sub-mesh for all incomplete elements
700  inline LagrangianSubMesh subIncomplete() const;
701 
702  //- Return the global sizes of all the sub-meshes. A value of -1
703  // will be set for processor-patch sub-meshes.
705 
706  //- Return the global positions
708 
709  //- Return the global position of an element
710  point position(const label i) const;
711 
712  //- Convert a position into a set of coordinates and a
713  // corresponding tetrahedron. Return a status flag indicating
714  // where the point is relative to the bounds of the mesh.
716  (
717  const point& position,
719  label& celli,
720  label& facei,
721  label& faceTrii,
722  const scalar fraction
723  ) const;
724 
725  //- Convert set of positions into a set of coordinates and a
726  // corresponding tetrahedron. Return status flags indicating
727  // where the points are relative to the bounds of the mesh.
729  (
730  const List<point>& position,
732  labelList& celli,
733  labelList& facei,
735  const scalarList& fraction
736  ) const;
737 
738 
739  // Modify
740 
741  //- Partition the mesh such that the groups are contiguous in memory
742  void partition();
743 
744  //- Track the positions along the given displacements
745  template<class Displacement>
746  void track
747  (
748  const List<LagrangianState>& endState,
749  const Displacement& displacement,
750  const LagrangianSubScalarField& deltaFraction,
752  );
753 
754  //- Cross the faces
755  void crossFaces
756  (
758  );
759 
760  //- Inject specified elements into the mesh. This method does not
761  // know how to set values for the fields. If any fields are
762  // registered then this function will error.
763  template<class ... FieldNamesAndFields>
765  (
767  const labelField& celli,
768  const labelField& facei,
769  const labelField& faceTrii,
770  const FieldNamesAndFields& ... fieldNamesAndFields
771  );
772 
773  //- Inject specified elements into the mesh. Fields are expanded
774  // accordingly. New field values can be specified manually by
775  // providing alternating name and field arguments. Source
776  // conditions will be used to set values for all other fields.
777  // Internal fields without source conditions must be specified
778  // manually or an error will result.
779  template<class ... FieldNamesAndFields>
781  (
782  const LagrangianInjection& injection,
784  const labelField& celli,
785  const labelField& facei,
786  const labelField& faceTrii,
787  const FieldNamesAndFields& ... fieldNamesAndFields
788  );
789 
790  //- Birth specified elements into the mesh. Fields are expanded
791  // accordingly. New field values can be specified manually by
792  // providing alternating name and field arguments. Values for
793  // other fields will be mapped from the parent elements.
794  template<class ... FieldNamesAndFields>
796  (
797  const labelList& parents,
798  const FieldNamesAndFields& ... fieldNamesAndFields
799  );
800 
801  //- Reset the mesh to the old-time conditions
802  void reset(const bool initial, const bool final);
803 
804  //- Clear all geometry out of the Lagrangian mesh
805  void clear();
806 
807 
808  // Advanced Modify
809 
810  //- Manually partition a sub-set of the mesh given a list of
811  // element-group pairs. Elements not given a group will remain in
812  // a block at the start of the mesh. Elements with multiple groups
813  // will be placed in the block associated with the group with the
814  // largest index. Return offsets to the blocks.
816  (
817  const label nGroups,
818  const UList<labelPair>& elementsGroups
819  );
820 
821  //- As above but for enumerated groups
822  template<class Enumeration>
824  (
825  const label nGroups,
826  const UList<elementGroup<Enumeration>>& elementsGroups
827  );
828 
829  //- Remove specified elements from the mesh. Shuffles everything
830  // else up.
831  void remove(const UList<label>& elements);
832 
833  //- Remove a specified number of elements from the end of the mesh.
834  void remove(const label nElements);
835 
836 
837  // Mesh changes
838 
839  //- Clear the positions uses during mapping
840  virtual void clearPosition();
841 
842  //- Store the positions for use during mapping
843  virtual void storePosition();
844 
845  //- Update topology using the given map
846  virtual void topoChange(const polyTopoChangeMap&);
847 
848  //- Update from another mesh using the given map
849  virtual void mapMesh(const polyMeshMap&);
850 
851  //- Redistribute or update using the given distribution map
852  virtual void distribute(const polyDistributionMap&);
853 
854 
855  // Write
856 
857  //- Write using given format, version and compression
858  virtual bool writeObject
859  (
863  const bool write = true
864  ) const;
865 
866  //- Write settings from the database
867  virtual bool write(const bool write = true) const;
868 
869 
870  // Member Operators
871 
872  //- Disallow default bitwise assignment
873  void operator=(const LagrangianMesh&) = delete;
874 };
875 
876 
877 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
878 
879 } // End namespace Foam
880 
881 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
882 
883 #include "LagrangianMeshI.H"
884 
885 #ifdef NoRepository
886  #include "LagrangianMeshTemplates.C"
887 #endif
888 
889 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
890 
891 #endif
892 
893 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Pre-declare SubField and related Field type.
Definition: Field.H:83
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
A HashTable with keys but without contents.
Definition: HashSet.H:62
An STL-conforming hash table.
Definition: HashTable.H:127
friend Ostream & operator(Ostream &, const HashTable< regIOobject *, word, string::hash > &)
readOption
Enumeration defining the read options.
Definition: IOobject.H:117
word group() const
Return group (extension part of name)
Definition: IOobject.C:321
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
Boundary part of a Lagrangian mesh. Just a list of Lagrangian patches with some added convenience fun...
Base class for Lagrangian source conditions.
Base class for Lagrangian injections. Minimal wrapper over LagrangianSource. Implements some utility ...
Class to define the scope of Lagrangian mesh state changes.
changer(LagrangianMesh &mesh, const LagrangianState state)
Construct for a Lagrangian mesh with a given state.
Class to hold element-group indices, and associate the group.
Class to hold and index into the field references associated with a.
const vector & operator()(const label i) const
The displacement for a given index.
linearDisplacement(const LagrangianSubVectorField &linear)
Construct from a reference to the displacement.
Class to hold and index into the field references associated with a.
Pair< vector > operator()(const label i) const
The displacements for a given index.
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.
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 LagrangianSubMesh & subAll() const
Return a sub-mesh for all elements.
LagrangianSubMesh inject(const barycentricField &coordinates, const labelField &celli, const labelField &facei, const labelField &faceTrii, const FieldNamesAndFields &... fieldNamesAndFields)
Inject specified elements into the mesh. This method does not.
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.
bool changing() const
Return whether or not the mesh is changing.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
static const NamedEnum< permutationAlgorithm, 2 > permutationAlgorithmNames_
Permutation algorithm names.
LagrangianBoundaryMesh BoundaryMesh
Boundary mesh type.
~LagrangianMesh()
Destructor.
LagrangianSubMesh birth(const labelList &parents, const FieldNamesAndFields &... fieldNamesAndFields)
Birth specified elements into the mesh. Fields are expanded.
ClassName("LagrangianMesh")
Run-time type information.
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.
LagrangianSubMesh subNone() const
Return a sub-mesh for no elements.
LagrangianSubMesh subIncomplete() const
Return a sub-mesh for all incomplete elements.
HashTable< const GeoField * > lookupCurrentFields(const bool strict=false) const
Lookup all current-time fields of the given type.
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.
label nStates() const
Return the number of states.
label subMeshIndex() const
Get a sub-mesh index.
LagrangianMesh Mesh
Mesh type.
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 operator=(const LagrangianMesh &)=delete
Disallow default bitwise assignment.
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.
label globalSize() const
Get the global number of elements.
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.
Base class for Lagrangian boundary conditions.
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...
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
const label & second() const
Return second.
Definition: PairI.H:121
const label & first() const
Return first.
Definition: PairI.H:107
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Centred interpolation interpolation scheme class.
Definition: linear.H:53
A non-conformal cyclic boundary condition for Lagrangian. Properties are transformed by the transform...
Registry of regIOobjects.
const Time & time() const
Return time.
const objectRegistry & thisDb() const
Return the object registry.
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
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A processor boundary condition for Lagrangian. Properties are communicated to and from the neighbour ...
Processor Lagrangian patch. This is used for the patches that interface between processors across fac...
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
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
LagrangianState
Lagrangian state enumeration.
LagrangianGroup
Lagrangian group enumeration.
labelList f(nPoints)