vtkPVblockMeshConvert.C 7.95 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

\*---------------------------------------------------------------------------*/

26
27
#include "vtkPVblockMesh.H"
#include "vtkPVblockMeshReader.h"
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45

// OpenFOAM includes
#include "blockMesh.H"
#include "Time.H"

#include "vtkOpenFOAMPoints.H"

// VTK includes
#include "vtkCellArray.h"
#include "vtkDataArraySelection.h"
#include "vtkMultiBlockDataSet.h"
#include "vtkPoints.h"
#include "vtkPolyData.h"
#include "vtkUnstructuredGrid.h"


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

46
void Foam::vtkPVblockMesh::convertMeshBlocks
47
48
49
50
51
52
53
54
55
56
57
(
    vtkMultiBlockDataSet* output,
    int& blockNo
)
{
    vtkDataArraySelection* selection = reader_->GetBlockSelection();
    arrayRange& range = arrayRangeBlocks_;
    range.block(blockNo);   // set output block
    label datasetNo = 0;       // restart at dataset 0

    const blockMesh& blkMesh = *meshPtr_;
58
    const Foam::pointField& blockPoints = blkMesh.vertices();
59
60
61

    if (debug)
    {
62
        Info<< "<beg> Foam::vtkPVblockMesh::convertMeshBlocks" << endl;
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
    }

    int blockI = 0;
    const scalar scaleFactor = blkMesh.scaleFactor();

    for
    (
        int partId = range.start();
        partId < range.end();
        ++partId, ++blockI
    )
    {
        if (!blockStatus_[partId])
        {
            continue;
        }

80
        const blockDescriptor& blockDef = blkMesh[blockI];
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132

        vtkUnstructuredGrid* vtkmesh = vtkUnstructuredGrid::New();

        // Convert OpenFOAM mesh vertices to VTK
        vtkPoints *vtkpoints = vtkPoints::New();
        vtkpoints->Allocate( blockDef.nPoints() );
        const labelList& blockLabels = blockDef.blockShape();

        vtkmesh->Allocate(1);
        vtkIdType nodeIds[8];

        forAll(blockLabels, ptI)
        {
            vtkInsertNextOpenFOAMPoint
            (
                vtkpoints,
                blockPoints[blockLabels[ptI]],
                scaleFactor
            );

            nodeIds[ptI] = ptI;
        }

        vtkmesh->InsertNextCell
        (
            VTK_HEXAHEDRON,
            8,
            nodeIds
        );

        vtkmesh->SetPoints(vtkpoints);
        vtkpoints->Delete();

        AddToBlock
        (
            output, vtkmesh, range, datasetNo,
            selection->GetArrayName(partId)
        );

        vtkmesh->Delete();
        datasetNo++;
    }


    // anything added?
    if (datasetNo)
    {
        ++blockNo;
    }

    if (debug)
    {
133
        Info<< "<end> Foam::vtkPVblockMesh::convertMeshBlocks" << endl;
134
135
136
137
    }
}


138
void Foam::vtkPVblockMesh::convertMeshEdges
139
140
141
142
143
144
145
146
147
148
149
150
(
    vtkMultiBlockDataSet* output,
    int& blockNo
)
{
    vtkDataArraySelection* selection = reader_->GetCurvedEdgesSelection();
    arrayRange& range = arrayRangeEdges_;

    range.block(blockNo);      // set output block
    label datasetNo = 0;       // restart at dataset 0

    const blockMesh& blkMesh = *meshPtr_;
151
    const blockEdgeList& edges = blkMesh.edges();
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170

    int edgeI = 0;
    const scalar scaleFactor = blkMesh.scaleFactor();

    for
    (
        int partId = range.start();
        partId < range.end();
        ++partId, ++edgeI
    )
    {
        if (!edgeStatus_[partId])
        {
            continue;
        }

        // search each block
        forAll(blkMesh, blockI)
        {
171
            const blockDescriptor& blockDef = blkMesh[blockI];
172
173
174

            edgeList blkEdges = blockDef.blockShape().edges();

175
            // List of edge point and weighting factors
176
            pointField edgesPoints[12];
177
178
179
            scalarList edgesWeights[12];
            blockDef.edgesPointsWeights(edgesPoints, edgesWeights);

180
181
182
183
184
185
186
187
188
189
190
191
192
            // find the corresponding edge within the block
            label foundEdgeI = -1;
            forAll(blkEdges, blkEdgeI)
            {
                if (edges[edgeI].compare(blkEdges[blkEdgeI]))
                {
                    foundEdgeI = blkEdgeI;
                    break;
                }
            }

            if (foundEdgeI != -1)
            {
193
                const List<point>& edgePoints = edgesPoints[foundEdgeI];
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245

                vtkPolyData* vtkmesh = vtkPolyData::New();
                vtkPoints* vtkpoints = vtkPoints::New();

                vtkpoints->Allocate( edgePoints.size() );
                vtkmesh->Allocate(1);

                vtkIdType pointIds[edgePoints.size()];
                forAll(edgePoints, ptI)
                {
                    vtkInsertNextOpenFOAMPoint
                    (
                        vtkpoints,
                        edgePoints[ptI],
                        scaleFactor
                    );
                    pointIds[ptI] = ptI;
                }

                vtkmesh->InsertNextCell
                (
                    VTK_POLY_LINE,
                    edgePoints.size(),
                    pointIds
                );

                vtkmesh->SetPoints(vtkpoints);
                vtkpoints->Delete();

                AddToBlock
                (
                    output, vtkmesh, range, datasetNo,
                    selection->GetArrayName(partId)
                );

                vtkmesh->Delete();
                datasetNo++;

                break;
            }
        }
    }


    // anything added?
    if (datasetNo)
    {
        ++blockNo;
    }

    if (debug)
    {
246
        Info<< "<end> Foam::vtkPVblockMesh::convertMeshEdges" << endl;
247
248
249
250
251
    }

}


252
void Foam::vtkPVblockMesh::convertMeshCorners
253
254
255
256
257
258
259
260
261
(
    vtkMultiBlockDataSet* output,
    int& blockNo
)
{
    arrayRange& range = arrayRangeCorners_;
    range.block(blockNo);      // set output block
    label datasetNo = 0;       // restart at dataset 0

262
    const pointField& blockPoints = meshPtr_->vertices();
263
264
265
266
    const scalar& scaleFactor = meshPtr_->scaleFactor();

    if (debug)
    {
267
        Info<< "<beg> Foam::vtkPVblockMesh::convertMeshCorners" << endl;
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
    }

    if (true)  // or some flag or other condition
    {
        vtkPolyData* vtkmesh = vtkPolyData::New();
        vtkPoints* vtkpoints = vtkPoints::New();
        vtkCellArray* vtkcells = vtkCellArray::New();

        vtkpoints->Allocate( blockPoints.size() );
        vtkcells->Allocate( blockPoints.size() );

        vtkIdType pointId = 0;
        forAll(blockPoints, ptI)
        {
            vtkInsertNextOpenFOAMPoint
            (
                vtkpoints,
                blockPoints[ptI],
                scaleFactor
            );

            vtkcells->InsertNextCell(1, &pointId);
            pointId++;
        }

        vtkmesh->SetPoints(vtkpoints);
        vtkpoints->Delete();

        vtkmesh->SetVerts(vtkcells);
        vtkcells->Delete();

        AddToBlock
        (
            output, vtkmesh, range, datasetNo,
            arrayRangeCorners_.name()
        );
        vtkmesh->Delete();

        datasetNo++;
    }

    // anything added?
    if (datasetNo)
    {
        ++blockNo;
    }

    if (debug)
    {
317
        Info<< "<end> Foam::vtkPVblockMesh::convertMeshCorners" << endl;
318
319
320
321
322
    }
}


// ************************************************************************* //