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-2018 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 "Cloud.H"
33 
34 #include "vtkOpenFOAMTupleRemap.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 template<class Type>
39 void Foam::vtkPVFoam::convertLagrangianFields
40 (
41  const IOobjectList& objects,
42  vtkMultiBlockDataSet* output,
43  const label datasetNo
44 )
45 {
46  const arrayRange& range = arrayRangeLagrangian_;
47 
48  forAllConstIter(IOobjectList, objects, iter)
49  {
50  // restrict to this IOField<Type>
51  if (iter()->headerClassName() == IOField<Type>::typeName)
52  {
53  IOField<Type> tf(*iter());
54  convertLagrangianField(tf, output, range, datasetNo);
55  }
56  }
57 }
58 
59 
60 template<class Type>
61 void Foam::vtkPVFoam::convertLagrangianField
62 (
63  const IOField<Type>& tf,
64  vtkMultiBlockDataSet* output,
65  const arrayRange& range,
66  const label datasetNo
67 )
68 {
69  const label nComp = pTraits<Type>::nComponents;
70 
71  vtkFloatArray* pointData = vtkFloatArray::New();
72  pointData->SetNumberOfTuples(tf.size());
73  pointData->SetNumberOfComponents(nComp);
74  pointData->Allocate(nComp*tf.size());
75  pointData->SetName(tf.name().c_str());
76 
77  if (debug)
78  {
79  Info<< "convert LagrangianField: "
80  << tf.name()
81  << " size = " << tf.size()
82  << " nComp=" << nComp
83  << " nTuples = " << tf.size() << endl;
84  }
85 
86  float vec[nComp];
87  forAll(tf, i)
88  {
89  const Type& t = tf[i];
90  for (direction d=0; d<nComp; ++d)
91  {
92  vec[d] = component(t, d);
93  }
94  vtkOpenFOAMTupleRemap<Type>(vec);
95 
96  pointData->InsertTuple(i, vec);
97  }
98 
99 
100  vtkPolyData::SafeDownCast
101  (
102  GetDataSetFromBlock(output, range, datasetNo)
103  ) ->GetPointData()
104  ->AddArray(pointData);
105 
106  pointData->Delete();
107 }
108 
109 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
110 
111 #endif
112 
113 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
static const char *const typeName
Definition: Field.H:105
uint8_t direction
Definition: direction.H:45
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const tensorField & tf
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
messageStream Info
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)