vtkPVFoamFields.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-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 "vtkPVFoam.H"
27 #include "vtkPVFoamReader.h"
28 #include "vtkPVFoamVolFields.H"
29 #include "vtkPVFoamPointFields.H"
31 
32 // OpenFOAM includes
33 #include "domainDecomposition.H"
34 #include "cloud.H"
35 #include "LagrangianMesh.H"
36 #include "IOobjectList.H"
37 
38 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39 
40 Foam::IOobjectList Foam::vtkPVFoam::getObjects
41 (
42  const wordHashSet& selected,
43  const fvMesh& mesh,
44  const fileName& instance,
45  const fileName& local
46 )
47 {
48  // If nothing is selected then return an empty list
49  if (selected.empty())
50  {
51  return IOobjectList(0);
52  }
53 
54  // Create the list of objects at the instance
55  IOobjectList objects(mesh, instance, local);
56 
57  // Add any objects from constant that are not already present
58  IOobjectList objectsConstant(mesh, Time::constant(), local);
59  forAllIter(IOobjectList, objectsConstant, iter)
60  {
61  if (!objects.found(iter.key()))
62  {
63  objects.add
64  (
65  *objectsConstant.HashPtrTable<IOobject>::remove(iter)
66  );
67  }
68  }
69 
70  // Remove everything that is not selected
71  forAllIter(IOobjectList, objects, iter)
72  {
73  if (!selected.found(iter()->name()))
74  {
75  objects.erase(iter);
76  }
77  }
78 
79  return objects;
80 }
81 
82 
83 void Foam::vtkPVFoam::convertFields(vtkMultiBlockDataSet* output)
84 {
85  if (reader_->GetDecomposedCase() && !procDbsPtr_->nProcs()) return;
86 
87  const wordHashSet selectedFields =
88  getSelected(reader_->GetFieldSelection());
89 
90  if (selectedFields.empty()) return;
91 
92  // Get objects (fields) for this time - only keep selected fields
93  // the region name is already in the mesh db
94  IOobjectList objects
95  (
96  getObjects
97  (
98  selectedFields,
99  reader_->GetDecomposedCase()
100  ? procMeshesPtr_->procMeshes().first()
101  : procMeshesPtr_->completeMesh(),
102  reader_->GetDecomposedCase()
103  ? procDbsPtr_->proc0Time().name()
104  : procDbsPtr_->completeTime().name()
105  )
106  );
107 
108  if (objects.empty()) return;
109 
111  << "Converting OpenFOAM volume fields" << endl;
112  forAllConstIter(IOobjectList, objects, iter)
113  {
114  DebugInfo
115  << " " << iter()->name()
116  << " == " << iter()->objectPath(false) << endl;
117  }
118 
119  PtrList<PrimitivePatchInterpolation<primitivePatch>>
120  ppInterpList(procMeshesPtr_->completeMesh().boundaryMesh().size());
121  forAll(ppInterpList, i)
122  {
123  ppInterpList.set
124  (
125  i,
126  new PrimitivePatchInterpolation<primitivePatch>
127  (
128  procMeshesPtr_->completeMesh().boundaryMesh()[i]
129  )
130  );
131  }
132 
133  bool interpFields = reader_->GetInterpolateVolFields();
134 
135  #define CONVERT_VOL_FIELDS(Type, nullArg) \
136  convertVolFields<Type> \
137  ( \
138  ppInterpList, \
139  objects, \
140  interpFields, \
141  output \
142  );
143  FOR_ALL_FIELD_TYPES(CONVERT_VOL_FIELDS)
144  #undef CONVERT_VOL_FIELDS
145 
146  #define CONVERT_VOL_INTERNAL_FIELDS(Type, nullArg) \
147  convertVolInternalFields<Type>(objects, output);
148  FOR_ALL_FIELD_TYPES(CONVERT_VOL_INTERNAL_FIELDS)
149  #undef CONVERT_VOL_INTERNAL_FIELDS
150 
151  #define CONVERT_SURFACE_FIELDS(Type, nullArg) \
152  convertSurfaceFields<Type>(objects, output);
153  FOR_ALL_FIELD_TYPES(CONVERT_SURFACE_FIELDS)
154  #undef CONVERT_SURFACE_FIELDS
155 
156  #define CONVERT_POINT_FIELDS(Type, nullArg) \
157  convertPointFields<Type>(objects, output);
158  FOR_ALL_FIELD_TYPES(CONVERT_POINT_FIELDS)
159  #undef CONVERT_POINT_FIELDS
160 }
161 
162 
163 void Foam::vtkPVFoam::convertlagrangianFields(vtkMultiBlockDataSet* output)
164 {
165  const fileName lagrangianPrefix =
166  meshRegion_ == polyMesh::defaultRegion
167  ? fileName(lagrangian::cloud::prefix)
168  : meshRegion_/lagrangian::cloud::prefix;
169 
170  arrayRange& range = arrayRangelagrangian_;
171 
172  const wordHashSet selectedFields = getSelected
173  (
174  reader_->GetlagrangianFieldSelection()
175  );
176 
177  if (selectedFields.empty()) return;
178 
179  for (int partId = range.start(); partId < range.end(); ++partId)
180  {
181  const word lagrangianName = getPartName(partId);
182  const label datasetNo = partDataset_[partId];
183 
184  if (!partStatus_[partId] || datasetNo < 0) continue;
185 
186  // Find a processor that has some of the cloud in it
187  label proci = -1;
188  if (reader_->GetDecomposedCase())
189  {
190  forAll(procDbsPtr_->procTimes(), procj)
191  {
192  if
193  (
195  (
196  fileHandler().filePath
197  (
198  procDbsPtr_->procTimes()[procj].path()
199  /procDbsPtr_().completeTime().name()
200  /lagrangianPrefix
201  /lagrangianName
202  )
203  )
204  )
205  {
206  proci = procj;
207  break;
208  }
209  }
210  }
211 
212  if (reader_->GetDecomposedCase() && proci == -1) continue;
213 
214  // Get the lagrangian fields for this time and this cloud
215  // but only keep selected fields
216  // the region name is already in the mesh db
217  IOobjectList objects
218  (
219  getObjects
220  (
221  selectedFields,
222  reader_->GetDecomposedCase()
223  ? procMeshesPtr_->procMeshes()[proci]
224  : procMeshesPtr_->completeMesh(),
225  procDbsPtr_().completeTime().name(),
226  lagrangian::cloud::prefix/lagrangianName
227  )
228  );
229 
230  if (objects.empty()) continue;
231 
233  << "Converting OpenFOAM lagrangian fields" << endl;
234  forAllConstIter(IOobjectList, objects, iter)
235  {
236  DebugInfo
237  << " " << iter()->name()
238  << " == " << iter()->objectPath(false) << endl;
239  }
240 
241  #define CONVERT_LAGRANGIAN_FIELDS(Type, nullArg) \
242  convertlagrangianFields<Type>(objects, output, datasetNo);
243  CONVERT_LAGRANGIAN_FIELDS(label, );
244  FOR_ALL_FIELD_TYPES(CONVERT_LAGRANGIAN_FIELDS);
245  #undef CONVERT_LAGRANGIAN_FIELDS
246  }
247 }
248 
249 
250 void Foam::vtkPVFoam::convertLagrangianFields(vtkMultiBlockDataSet* output)
251 {
252  arrayRange& range = arrayRangeLagrangian_;
253 
254  const wordHashSet selectedFields = getSelected
255  (
256  reader_->GetLagrangianFieldSelection()
257  );
258 
259  if (selectedFields.empty()) return;
260 
261  for (int partId = range.start(); partId < range.end(); ++partId)
262  {
263  const word LagrangianName = getPartName(partId);
264  const label datasetNo = partDataset_[partId];
265 
266  if (!partStatus_[partId] || datasetNo < 0) continue;
267 
268  // Get the Lagrangian fields for this time and this cloud
269  // but only keep selected fields
270  // the region name is already in the mesh db
271  IOobjectList objects
272  (
273  getObjects
274  (
275  selectedFields,
276  reader_->GetDecomposedCase()
277  ? procMeshesPtr_->procMeshes().first()
278  : procMeshesPtr_->completeMesh(),
279  procDbsPtr_().completeTime().name(),
280  LagrangianMesh::prefix/LagrangianName
281  )
282  );
283 
284  if (objects.empty()) continue;
285 
287  << "Converting OpenFOAM Lagrangian fields" << endl;
288  forAllConstIter(IOobjectList, objects, iter)
289  {
290  DebugInfo
291  << " " << iter()->name()
292  << " == " << iter()->objectPath(false) << endl;
293  }
294 
295  #define CONVERT_LAGRANGIAN_FIELDS(Type, GeoField) \
296  convertLagrangianFields<Type, GeoField>(objects, output, datasetNo);
297  CONVERT_LAGRANGIAN_FIELDS(label, LagrangianField);
298  FOR_ALL_FIELD_TYPES(CONVERT_LAGRANGIAN_FIELDS, LagrangianField);
299  CONVERT_LAGRANGIAN_FIELDS(label, LagrangianInternalField);
300  FOR_ALL_FIELD_TYPES(CONVERT_LAGRANGIAN_FIELDS, LagrangianInternalField);
301  #undef CONVERT_LAGRANGIAN_FIELDS
302  }
303 }
304 
305 
306 // ************************************************************************* //
scalar range
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:458
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
static const word prefix
Instance prefix.
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:69
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:270
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
const fileOperation & fileHandler()
Get current file handler.
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
DimensionedField< Type, LagrangianMesh > LagrangianInternalField
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
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
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
GeometricField< Type, LagrangianMesh > LagrangianField
objects