initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / applications / utilities / postProcessing / graphics / PVFoamReader / vtkFoam / vtkFoamConvertPointField.H
blobe199e10f0321ec02f9860ffb66bb442844743b16
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 InClass
26     Foam::vtkFoam
28 \*---------------------------------------------------------------------------*/
30 #ifndef vtkFoamConvertPointField_H
31 #define vtkFoamConvertPointField_H
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 #include "interpolatePointToCell.H"
37 template<class Type>
38 void Foam::vtkFoam::convertPointField
40     const GeometricField<Type, pointPatchField, pointMesh>& ptf,
41     const GeometricField<Type, fvPatchField, volMesh>& tf
44     vtkUnstructuredGrid *vtkMesh =
45         vtkUnstructuredGrid::SafeDownCast(reader_->GetOutput(0));
47     vtkFloatArray *pointTypes = vtkFloatArray::New();
48     pointTypes->SetNumberOfTuples(ptf.size() + addPointCellLabels_.size());
49     pointTypes->SetNumberOfComponents(Type::nComponents);
50     pointTypes->Allocate(Type::nComponents*ptf.size());
51     pointTypes->SetName(ptf.name().c_str());
53     float vec[Type::nComponents];
55     forAll(ptf, i)
56     {
57         for (direction d=0; d<Type::nComponents; d++)
58         {
59             vec[d] = ptf[i][d];
60         }
62         pointTypes->InsertTuple(i, vec);
63     }
65     label i = ptf.size();
67     if (&tf != &GeometricField<Type, fvPatchField, volMesh>::null())
68     {
69         forAll(addPointCellLabels_, api)
70         {
71             Type t = tf[addPointCellLabels_[api]];
73             for (direction d=0; d<Type::nComponents; d++)
74             {
75                 vec[d] = t[d];
76             }
78             pointTypes->InsertTuple(i++, vec);
79         }
80     }
81     else
82     {
83         forAll(addPointCellLabels_, api)
84         {
85             Type t = interpolatePointToCell(ptf, addPointCellLabels_[api]);
87             for (direction d=0; d<Type::nComponents; d++)
88             {
89                 vec[d] = t[d];
90             }
92             pointTypes->InsertTuple(i++, vec);
93         }
94     }
96     vtkMesh->GetPointData()->AddArray(pointTypes);
97     pointTypes->Delete();
101 template<>
102 void Foam::vtkFoam::convertPointField
104     const GeometricField<scalar, pointPatchField, pointMesh>& psf,
105     const GeometricField<scalar, fvPatchField, volMesh>& sf
108     vtkUnstructuredGrid *vtkMesh =
109         vtkUnstructuredGrid::SafeDownCast(reader_->GetOutput(0));
111     vtkFloatArray *pointScalars = vtkFloatArray::New();
112     pointScalars->SetNumberOfTuples(psf.size() + addPointCellLabels_.size());
113     pointScalars->SetNumberOfComponents(1);
114     pointScalars->Allocate(psf.size());
115     pointScalars->SetName(psf.name().c_str());
117     for (int i=0; i<psf.size(); i++)
118     {
119         pointScalars->InsertComponent(i, 0, psf[i]);
120     }
122     label i = psf.size();
124     if (&sf != &GeometricField<scalar, fvPatchField, volMesh>::null())
125     {
126         forAll(addPointCellLabels_, api)
127         {
128             pointScalars->InsertComponent
129             (
130                 i++,
131                 0,
132                 sf[addPointCellLabels_[api]]
133             );
134         }
135     }
136     else
137     {
138         forAll(addPointCellLabels_, api)
139         {
140             pointScalars->InsertComponent
141             (
142                 i++,
143                 0,
144                 interpolatePointToCell(psf, addPointCellLabels_[api])
145             );
146         }
147     }
149     vtkMesh->GetPointData()->AddArray(pointScalars);
150     if (!vtkMesh->GetPointData()->GetScalars())
151     {
152         vtkMesh->GetPointData()->SetScalars(pointScalars);
153     }
155     pointScalars->Delete();
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161 #endif
163 // ************************************************************************* //