vtkPVFoamLagrangianFields.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-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 InClass
25  vtkPVFoam
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #ifndef vtkPVFoamLagrangianFields_H
30 #define vtkPVFoamLagrangianFields_H
31 
32 #include "vtkPVFoam.H"
33 #include "vtkPVFoamReader.h"
34 #include "vtkOpenFOAMTupleRemap.H"
35 
36 // OpenFOAM includes
37 #include "LagrangianFields.H"
40 
41 // VTK includes
42 #include "vtkFloatArray.h"
43 #include "vtkPointData.h"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 template<class Type>
48 void Foam::vtkPVFoam::convertlagrangianFields
49 (
50  const IOobjectList& objects,
51  vtkMultiBlockDataSet* output,
52  const label datasetNo
53 )
54 {
55  const arrayRange& range = arrayRangelagrangian_;
56 
57  forAllConstIter(IOobjectList, objects, iter)
58  {
59  // restrict to this IOField<Type>
60  if (iter()->headerClassName() != IOField<Type>::typeName) continue;
61 
62  const IOField<Type> tf
63  (
64  reader_->GetDecomposedCase()
65  ? lagrangianReconstructors_[datasetNo]
66  .reconstructField<Type, IOField, IOField>(*iter())
67  : tmp<IOField<Type>>(new IOField<Type>(*iter()))
68  );
69 
70  convertlagrangianField
71  (
72  tf,
73  output,
74  range,
75  datasetNo
76  );
77  }
78 }
79 
80 
81 template<class Type, template<class> class GeoField>
82 void Foam::vtkPVFoam::convertLagrangianFields
83 (
84  const IOobjectList& objects,
85  vtkMultiBlockDataSet* output,
86  const label datasetNo
87 )
88 {
89  const arrayRange& range = arrayRangeLagrangian_;
90 
91  forAllConstIter(IOobjectList, objects, iter)
92  {
93  // restrict to this GeoField<Type>
94  if (iter()->headerClassName() != GeoField<Type>::typeName) continue;
95 
96  const GeoField<Type> tf
97  (
98  reader_->GetDecomposedCase()
99  ? LagrangianReconstructors_[datasetNo]
100  .reconstructField<GeoField<Type>>(*iter())
101  : tmp<GeoField<Type>>
102  (
103  new GeoField<Type>(*iter(), LagrangianMeshes_[datasetNo])
104  )
105  );
106 
107  convertLagrangianField<Type, GeoField>
108  (
109  tf,
110  output,
111  range,
112  datasetNo
113  );
114  }
115 }
116 
117 
118 template<class Type>
119 void Foam::vtkPVFoam::convertlagrangianField
120 (
121  const IOField<Type>& tf,
122  vtkMultiBlockDataSet* output,
123  const arrayRange& range,
124  const label datasetNo
125 )
126 {
127  const label nComp = pTraits<Type>::nComponents;
128 
129  vtkFloatArray* pointData = vtkFloatArray::New();
130  pointData->SetNumberOfTuples(tf.size());
131  pointData->SetNumberOfComponents(nComp);
132  pointData->Allocate(nComp*tf.size());
133  pointData->SetName(tf.name().c_str());
134 
136  << "Converting lagrangian field: " << tf.name()
137  << " size=" << tf.size()
138  << ", nComp=" << nComp << endl;
139 
140  float vec[nComp];
141 
142  forAll(tf, i)
143  {
144  const Type& t = tf[i];
145 
146  for (direction d=0; d<nComp; ++d)
147  {
148  vec[d] = component(t, d);
149  }
150 
151  vtkOpenFOAMTupleRemap<Type>(vec);
152 
153  pointData->InsertTuple(i, vec);
154  }
155 
156  vtkPolyData::SafeDownCast
157  (
158  GetDataSetFromBlock(output, range, datasetNo)
159  ) ->GetPointData()
160  ->AddArray(pointData);
161 
162  pointData->Delete();
163 }
164 
165 
166 template<class Type, template<class> class GeoField>
167 void Foam::vtkPVFoam::convertLagrangianField
168 (
169  const GeoField<Type>& tf,
170  vtkMultiBlockDataSet* output,
171  const arrayRange& range,
172  const label datasetNo
173 )
174 {
175  static const label nComp = pTraits<Type>::nComponents;
176 
177  vtkFloatArray* pointData = vtkFloatArray::New();
178  pointData->SetNumberOfTuples(tf.primitiveField().size());
179  pointData->SetNumberOfComponents(nComp);
180  pointData->Allocate(nComp*tf.primitiveField().size());
181  pointData->SetName(tf.name().c_str());
182 
184  << "Converting Lagrangian field: " << tf.name()
185  << " size=" << tf.primitiveField().size()
186  << ", nComp=" << nComp << endl;
187 
188  float vec[nComp];
189 
190  forAll(tf.primitiveField(), i)
191  {
192  const Type& t = tf.primitiveField()[i];
193 
194  for (direction d=0; d<nComp; ++d)
195  {
196  vec[d] = component(t, d);
197  }
198 
199  vtkOpenFOAMTupleRemap<Type>(vec);
200 
201  pointData->InsertTuple(i, vec);
202  }
203 
204  vtkPolyData::SafeDownCast
205  (
206  GetDataSetFromBlock(output, range, datasetNo)
207  ) ->GetPointData()
208  ->AddArray(pointData);
209 
210  pointData->Delete();
211 }
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 #endif
216 
217 // ************************************************************************* //
scalar range
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
static const char *const typeName
Definition: Field.H:106
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const tensorField & tf
#define DebugInFunction
Report an information message using Foam::Info.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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
void component(LagrangianPatchField< typename LagrangianPatchField< Type >::cmptType > &sf, const LagrangianPatchField< Type > &f, const direction d)
uint8_t direction
Definition: direction.H:45
objects