forked from ufz/ogs
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathVtkMappedMesh.cpp
183 lines (154 loc) · 4.65 KB
/
VtkMappedMesh.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
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
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
/**
* \file
* \author Lars Bilke
* \date 2014-02-27
* \brief Implementation of the VtkMappedMesh class.
*
* \copyright
* Copyright (c) 2012-2017, OpenGeoSys Community (http://www.opengeosys.org)
* Distributed under a Modified BSD License.
* See accompanying file LICENSE.txt or
* http://www.opengeosys.org/project/license
*
*/
#include "VtkMappedMesh.h"
#include <algorithm>
#include <logog/include/logog.hpp>
#include <vtkCellType.h>
#include <vtkCellTypes.h>
#include <vtkGenericCell.h>
#include <vtkIdTypeArray.h>
#include <vtkObjectFactory.h>
#include <vtkPoints.h>
#include "MeshLib/Elements/Element.h"
#include "MeshLib/MeshEnums.h"
#include "MeshLib/Node.h"
#include "MeshLib/VtkOGSEnum.h"
namespace MeshLib {
vtkStandardNewMacro(VtkMappedMesh)
vtkStandardNewMacro(VtkMappedMeshImpl)
void VtkMappedMeshImpl::PrintSelf(ostream &os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
os << indent << "NumberOfCells: " << this->NumberOfCells << endl;
}
void VtkMappedMeshImpl::SetNodes(std::vector<MeshLib::Node*> const & nodes)
{
_nodes = &nodes;
}
void VtkMappedMeshImpl::SetElements(std::vector<MeshLib::Element*> const & elements)
{
_elements = &elements;
this->NumberOfCells = _elements->size();
this->Modified();
}
vtkIdType VtkMappedMeshImpl::GetNumberOfCells()
{
return this->NumberOfCells;
}
int VtkMappedMeshImpl::GetCellType(vtkIdType cellId)
{
return OGSToVtkCellType((*_elements)[cellId]->getCellType());
}
void VtkMappedMeshImpl::GetCellPoints(vtkIdType cellId, vtkIdList *ptIds)
{
const MeshLib::Element* const elem = (*_elements)[cellId];
const unsigned numNodes(elem->getNumberOfNodes());
const MeshLib::Node* const* nodes = (*_elements)[cellId]->getNodes();
ptIds->SetNumberOfIds(numNodes);
for (unsigned i(0); i < numNodes; ++i)
ptIds->SetId(i, nodes[i]->getID());
if(GetCellType(cellId) == VTK_WEDGE)
{
for (unsigned i=0; i<3; ++i)
{
const unsigned prism_swap_id = ptIds->GetId(i);
ptIds->SetId(i, ptIds->GetId(i+3));
ptIds->SetId(i+3, prism_swap_id);
}
}
else if(GetCellType(cellId) == VTK_QUADRATIC_WEDGE)
{
std::array<vtkIdType, 15> ogs_nodeIds;
for (unsigned i=0; i<15; ++i)
ogs_nodeIds[i] = ptIds->GetId(i);
for (unsigned i=0; i<3; ++i)
{
ptIds->SetId(i, ogs_nodeIds[i+3]);
ptIds->SetId(i+3, ogs_nodeIds[i]);
}
for (unsigned i=0; i<3; ++i)
ptIds->SetId(6+i, ogs_nodeIds[8-i]);
for (unsigned i=0; i<3; ++i)
ptIds->SetId(9+i, ogs_nodeIds[14-i]);
ptIds->SetId(12, ogs_nodeIds[9]);
ptIds->SetId(13, ogs_nodeIds[11]);
ptIds->SetId(14, ogs_nodeIds[10]);
}
}
void VtkMappedMeshImpl::GetPointCells(vtkIdType ptId, vtkIdList *cellIds)
{
cellIds->Reset();
auto elements((*_nodes)[ptId]->getElements());
for (auto elem(elements.begin()); elem != elements.end(); ++elem)
cellIds->InsertNextId((*elem)->getID());
}
int VtkMappedMeshImpl::GetMaxCellSize()
{
unsigned int size = 0;
for (auto elem(_elements->begin()); elem != _elements->end(); ++elem)
size = std::max(size, (*elem)->getNumberOfNodes());
return size;
}
void VtkMappedMeshImpl::GetIdsOfCellsOfType(int type, vtkIdTypeArray *array)
{
array->Reset();
for (auto elem(_elements->begin()); elem != _elements->end(); ++elem)
{
if ((*elem)->getCellType() == VtkCellTypeToOGS(type))
array->InsertNextValue((*elem)->getID());
}
}
int VtkMappedMeshImpl::IsHomogeneous()
{
MeshLib::CellType type = (*(_elements->begin()))->getCellType();
for (auto elem(_elements->begin()); elem != _elements->end(); ++elem)
if((*elem)->getCellType() != type)
return 0;
return 1;
}
void VtkMappedMeshImpl::Allocate(vtkIdType, int)
{
vtkErrorMacro("Read only container.")
return;
}
vtkIdType VtkMappedMeshImpl::InsertNextCell(int, vtkIdList*)
{
vtkErrorMacro("Read only container.")
return -1;
}
vtkIdType VtkMappedMeshImpl::InsertNextCell(int, vtkIdType, vtkIdType*)
{
vtkErrorMacro("Read only container.")
return -1;
}
vtkIdType VtkMappedMeshImpl::InsertNextCell(
int, vtkIdType, vtkIdType*, vtkIdType, vtkIdType*)
{
vtkErrorMacro("Read only container.")
return -1;
}
void VtkMappedMeshImpl::ReplaceCell(vtkIdType, int, vtkIdType*)
{
vtkErrorMacro("Read only container.")
return;
}
VtkMappedMeshImpl::VtkMappedMeshImpl()
: _nodes(NULL), _elements(NULL), NumberOfCells(0)
{
}
VtkMappedMeshImpl::~VtkMappedMeshImpl()
{
// delete [] this->Elements;
}
} // end namespace