lagrangianFieldDecomposer.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 
55 class lagrangianFieldDecomposer
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
67 
68 
69  // Private Member Functions
70 
71  //- Disallow default bitwise copy construct
73 
74  //- Disallow default bitwise assignment
76 
77 
78 public:
79 
80  // Constructors
81 
82  //- Construct from components
84  (
85  const polyMesh& mesh,
86  const polyMesh& procMesh,
88  const labelList& cellProcAddressing,
89  const word& cloudName,
90  const Cloud<indexedParticle>& lagrangianPositions,
91  const List<SLList<indexedParticle*>*>& cellParticles
92  );
93 
94 
95  // Member Functions
96 
97  // Read the fields and hold on the pointer list
98  template<class Type>
99  static void readFields
100  (
101  const label cloudI,
102  const IOobjectList& lagrangianObjects,
103  PtrList<PtrList<IOField<Type> > >& lagrangianFields
104 // PtrList<IOField<Type> >& lagrangianFields
105  );
106 
107  template<class Type>
108  static void readFieldFields
109  (
110  const label cloudI,
111  const IOobjectList& lagrangianObjects,
112  PtrList
113  <
114  PtrList<CompactIOField<Field<Type>, Type> >
115  >& lagrangianFields
116 // PtrList<CompactIOField<Field<Type>, Type > >& lagrangianFields
117  );
118 
119 
120  //- Decompose volume field
121  template<class Type>
122  tmp<IOField<Type> > decomposeField
123  (
124  const word& cloudName,
125  const IOField<Type>& field
126  ) const;
127 
128  template<class Type>
129  tmp<CompactIOField<Field<Type>, Type> > decomposeFieldField
130  (
131  const word& cloudName,
132  const CompactIOField<Field<Type>, Type>& field
133  ) const;
134 
135 
136  template<class GeoField>
137  void decomposeFields
138  (
139  const word& cloudName,
140  const PtrList<GeoField>& fields
141  ) const;
142 
143  template<class GeoField>
145  (
146  const word& cloudName,
147  const PtrList<GeoField>& fields
148  ) const;
149 };
150 
151 
152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 
154 } // End namespace Foam
155 
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 
158 #ifdef NoRepository
159  #include "lagrangianFieldDecomposerDecomposeFields.C"
160 #endif
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 #endif
165 
166 // ************************************************************************* //
Foam::lagrangianFieldDecomposer::operator=
void operator=(const lagrangianFieldDecomposer &)
Disallow default bitwise assignment.
Foam::lagrangianFieldDecomposer::decomposeField
tmp< IOField< Type > > decomposeField(const word &cloudName, const IOField< Type > &field) const
Decompose volume field.
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::lagrangianFieldDecomposer::decomposeFieldField
tmp< CompactIOField< Field< Type >, Type > > decomposeFieldField(const word &cloudName, const CompactIOField< Field< Type >, Type > &field) const
fields
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p - rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:127
Foam::lagrangianFieldDecomposer::decomposeFieldFields
void decomposeFieldFields(const word &cloudName, const PtrList< GeoField > &fields) const
Definition: lagrangianFieldDecomposerDecomposeFields.C:203
Cloud.H
Foam::lagrangianFieldDecomposer::decomposeFields
void decomposeFields(const word &cloudName, const PtrList< GeoField > &fields) const
Definition: lagrangianFieldDecomposerDecomposeFields.C:186
Foam::lagrangianFieldDecomposer::particleIndices_
labelList particleIndices_
The indices of the particles on this processor.
Definition: lagrangianFieldDecomposer.H:65
Foam::label
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
passiveParticle.H
Foam::lagrangianFieldDecomposer::readFieldFields
static void readFieldFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< CompactIOField< Field< Type >, Type > > > &lagrangianFields)
Definition: lagrangianFieldDecomposerDecomposeFields.C:68
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
indexedParticle.H
Foam::lagrangianFieldDecomposer::positions_
Cloud< passiveParticle > positions_
Lagrangian positions for this processor.
Definition: lagrangianFieldDecomposer.H:62
Foam::lagrangianFieldDecomposer::readFields
static void readFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< IOField< Type > > > &lagrangianFields)
Definition: lagrangianFieldDecomposerDecomposeFields.C:33
Foam::lagrangianFieldDecomposer::procMesh_
const polyMesh & procMesh_
Reference to processor mesh.
Definition: lagrangianFieldDecomposer.H:59
Foam::lagrangianFieldDecomposer::lagrangianFieldDecomposer
lagrangianFieldDecomposer(const lagrangianFieldDecomposer &)
Disallow default bitwise copy construct.
cloudName
const word cloudName(propsDict.lookup("cloudName"))
faceProcAddressing
PtrList< labelIOList > & faceProcAddressing
Definition: checkFaceAddressingComp.H:9
List
Definition: Test.C:19
CompactIOField.H