lagrangianFieldDecomposer.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) 2011-2022 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::lagrangianFieldDecomposer
26 
27 Description
28  Lagrangian field decomposer.
29 
30 SourceFiles
31  lagrangianFieldDecomposer.C
32  lagrangianFieldDecomposerDecomposeFields.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef lagrangianFieldDecomposer_H
37 #define lagrangianFieldDecomposer_H
38 
39 #include "Cloud.H"
40 #include "CompactIOField.H"
41 #include "indexedParticle.H"
42 #include "passiveParticle.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 class IOobjectList;
50 
51 /*---------------------------------------------------------------------------*\
52  Class lagrangianFieldDecomposer Declaration
53 \*---------------------------------------------------------------------------*/
54 
56 {
57  // Private Data
58 
59  //- Reference to processor mesh
60  const polyMesh& procMesh_;
61 
62  //- Lagrangian positions for this processor
63  Cloud<passiveParticle> positions_;
64 
65  //- The indices of the particles on this processor
66  labelList particleIndices_;
67 
68 
69 public:
70 
71  // Constructors
72 
73  //- Construct from components
75  (
76  const polyMesh& mesh,
77  const polyMesh& procMesh,
78  const labelList& faceProcAddressing,
79  const labelList& cellProcAddressing,
80  const word& cloudName,
81  const Cloud<indexedParticle>& lagrangianPositions,
82  const List<SLList<indexedParticle*>*>& cellParticles
83  );
84 
85  //- Disallow default bitwise copy construction
87 
88 
89  // Member Functions
90 
91  // Read the fields and hold on the pointer list
92  template<class Type>
93  static void readFields
94  (
95  const label cloudI,
96  const IOobjectList& lagrangianObjects,
97  PtrList<PtrList<IOField<Type>>>& lagrangianFields
98  );
99 
100  template<class Type>
101  static void readFieldFields
102  (
103  const label cloudI,
104  const IOobjectList& lagrangianObjects,
105  PtrList<PtrList<CompactIOField<Field<Type>>>>& lagrangianFields
106  );
107 
108 
109  //- Decompose volume field
110  template<class Type>
112  (
113  const word& cloudName,
114  const IOField<Type>& field
115  ) const;
116 
117  template<class Type>
119  (
120  const word& cloudName,
121  const CompactIOField<Field<Type>>& field
122  ) const;
123 
124 
125  template<class GeoField>
126  void decomposeFields
127  (
128  const word& cloudName,
130  ) const;
131 
132  template<class GeoField>
134  (
135  const word& cloudName,
137  ) const;
138 
139 
140  // Member Operators
141 
142  //- Disallow default bitwise assignment
143  void operator=(const lagrangianFieldDecomposer&) = delete;
144 };
145 
146 
147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 
149 } // End namespace Foam
150 
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 
153 #ifdef NoRepository
155 #endif
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 #endif
160 
161 // ************************************************************************* //
A Field of objects of type <Type> with automated input and output using a compact storage....
Pre-declare SubField and related Field type.
Definition: Field.H:82
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
Template class for non-intrusive linked lists.
Definition: LList.H:76
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
void operator=(const lagrangianFieldDecomposer &)=delete
Disallow default bitwise assignment.
static void readFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< IOField< Type >>> &lagrangianFields)
tmp< IOField< Type > > decomposeField(const word &cloudName, const IOField< Type > &field) const
Decompose volume field.
lagrangianFieldDecomposer(const polyMesh &mesh, const polyMesh &procMesh, const labelList &faceProcAddressing, const labelList &cellProcAddressing, const word &cloudName, const Cloud< indexedParticle > &lagrangianPositions, const List< SLList< indexedParticle * > * > &cellParticles)
Construct from components.
static void readFieldFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< CompactIOField< Field< Type >>>> &lagrangianFields)
tmp< CompactIOField< Field< Type > > > decomposeFieldField(const word &cloudName, const CompactIOField< Field< Type >> &field) const
void decomposeFields(const word &cloudName, const PtrList< GeoField > &fields) const
void decomposeFieldFields(const word &cloudName, const PtrList< GeoField > &fields) const
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word cloudName(propsDict.lookup("cloudName"))