Skip to content
GitLab
Menu
Projects
Groups
Snippets
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
iwr
dune-vtk
Commits
3538348a
Commit
3538348a
authored
Aug 24, 2018
by
Praetorius, Simon
Browse files
support for polygonal elements
parent
eb544147
Changes
5
Hide whitespace changes
Inline
Side-by-side
dune.module
View file @
3538348a
...
...
@@ -8,3 +8,4 @@ Version: 0.1
Maintainer
:
simon
.
praetorius
@
tu
-
dresden
.
de
#depending on
Depends
:
dune
-
grid
dune
-
functions
Suggests
:
dune
-
uggrid
dune
-
polygongrid
dune/vtk/vtktypes.cc
View file @
3538348a
...
...
@@ -82,6 +82,14 @@ CellType::CellType (GeometryType const& t, CellParametrization parametrization)
permutation_
=
{
0
,
1
,
3
,
2
,
4
};
noPermutation_
=
false
;
}
else
if
(
t
.
isNone
()
&&
t
.
dim
()
==
1
)
{
type_
=
LINE
;
permutation_
=
{
0
,
1
};
}
else
if
(
t
.
isNone
()
&&
t
.
dim
()
==
2
)
{
type_
=
POLYGON
;
permutation_
=
{
0
,
1
,
2
,
3
,
4
,
5
,
6
,
7
,
8
,
9
,
10
,
11
,
12
,
13
,
14
,
15
,
16
,
17
,
18
,
19
};
}
else
{
std
::
cerr
<<
"Geometry Type not supported by VTK!
\n
"
;
std
::
abort
();
...
...
dune/vtk/vtktypes.hh
View file @
3538348a
...
...
@@ -36,16 +36,16 @@ namespace Dune { namespace experimental
enum
CellTypes
:
std
::
uint8_t
{
// Linear VTK cell types
VERTEX
=
1
,
POLY_VERTEX
=
2
,
// not supported
/*
POLY_VERTEX = 2, // not supported
*/
LINE
=
3
,
POLY_LINE
=
4
,
// not supported
/*
POLY_LINE = 4, // not supported
*/
TRIANGLE
=
5
,
TRIANGLE_STRIP
=
6
,
// not supported
POLYGON
=
7
,
// not supported
PIXEL
=
8
,
// not supported
/*
TRIANGLE_STRIP = 6, // not supported
*/
POLYGON
=
7
,
/*
PIXEL = 8, // not supported
*/
QUAD
=
9
,
TETRA
=
10
,
VOXEL
=
11
,
// not supported
/*
VOXEL = 11, // not supported
*/
HEXAHEDRON
=
12
,
WEDGE
=
13
,
PYRAMID
=
14
,
...
...
src/CMakeLists.txt
View file @
3538348a
...
...
@@ -18,4 +18,10 @@ add_executable("quadratic" quadratic.cc)
target_link_dune_default_libraries
(
"quadratic"
)
target_link_libraries
(
"quadratic"
dunevtk
)
add_subdirectory
(
test
)
\ No newline at end of file
if
(
dune-polygongrid_FOUND
)
add_executable
(
"polygongrid"
polygongrid.cc
)
target_link_dune_default_libraries
(
"polygongrid"
)
target_link_libraries
(
"polygongrid"
dunevtk dunepolygongrid
)
endif
()
add_subdirectory
(
test
)
src/polygongrid.cc
0 → 100644
View file @
3538348a
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifdef HAVE_CONFIG_H
# include "config.h"
#endif
#include <iostream>
#include <vector>
#include <dune/common/parallel/mpihelper.hh> // An initializer of MPI
#include <dune/common/exceptions.hh> // We use exceptions
#include <dune/common/filledarray.hh>
#include <dune/functions/functionspacebases/defaultglobalbasis.hh>
#include <dune/functions/functionspacebases/lagrangebasis.hh>
#include <dune/functions/functionspacebases/interpolate.hh>
#include <dune/functions/gridfunctions/analyticgridviewfunction.hh>
#include <dune/functions/gridfunctions/discreteglobalbasisfunction.hh>
#include <dune/polygongrid/grid.hh>
#include <dune/polygongrid/gridfactory.hh>
#include <dune/vtk/vtkwriter.hh>
using
namespace
Dune
;
using
namespace
Dune
::
experimental
;
using
namespace
Dune
::
Functions
;
using
GridType
=
Dune
::
PolygonGrid
<
double
>
;
std
::
unique_ptr
<
GridType
>
createArbitraryGrid
()
{
const
std
::
vector
<
Dune
::
FieldVector
<
double
,
2
>
>
vertices
=
{
{
0.0
,
0.0
},
{
0.5
,
0.0
},
{
1.0
,
0.0
},
{
0.0
,
0.4
},
{
0.5
,
0.2
},
{
0.7
,
0.4
},
{
1.0
,
0.4
},
{
0.0
,
0.7
},
{
0.7
,
0.6
},
{
1.0
,
0.6
},
{
0.0
,
1.0
},
{
0.3
,
1.0
},
{
1.0
,
1.0
}
};
const
std
::
vector
<
std
::
vector
<
unsigned
int
>
>
polys
=
{
{
0
,
1
,
4
,
3
},
{
1
,
2
,
6
,
5
,
4
},
{
3
,
4
,
5
,
8
,
11
,
7
},
{
5
,
6
,
9
,
8
},
{
7
,
11
,
10
},
{
8
,
9
,
12
,
11
}
};
Dune
::
GridFactory
<
GridType
>
factory
;
for
(
const
auto
&
vertex
:
vertices
)
factory
.
insertVertex
(
vertex
);
for
(
const
auto
&
poly
:
polys
)
factory
.
insertElement
(
Dune
::
GeometryTypes
::
none
(
2
),
poly
);
return
std
::
unique_ptr
<
GridType
>
(
factory
.
createGrid
()
);
}
int
main
(
int
argc
,
char
**
argv
)
{
Dune
::
MPIHelper
::
instance
(
argc
,
argv
);
auto
gridPtr
=
createArbitraryGrid
();
using
GridView
=
typename
GridType
::
LeafGridView
;
GridView
gridView
=
gridPtr
->
leafGridView
();
using
Writer
=
VtkWriter
<
GridView
>
;
Writer
vtkWriter
(
gridView
);
auto
p1Analytic
=
makeAnalyticGridViewFunction
([](
auto
const
&
x
)
{
return
std
::
sin
(
10
*
x
[
0
])
*
std
::
cos
(
10
*
x
[
1
])
+
std
::
sin
(
10
*
x
[
2
]);
},
gridView
);
vtkWriter
.
addPointData
(
p1Analytic
,
"p1"
);
vtkWriter
.
addCellData
(
p1Analytic
,
"p0"
);
vtkWriter
.
write
(
"poly_ascii_float32.vtu"
,
Vtk
::
ASCII
);
vtkWriter
.
write
(
"poly_binary_float32.vtu"
,
Vtk
::
BINARY
);
vtkWriter
.
write
(
"poly_compressed_float32.vtu"
,
Vtk
::
COMPRESSED
);
vtkWriter
.
write
(
"poly_ascii_float64.vtu"
,
Vtk
::
ASCII
,
Vtk
::
FLOAT64
);
vtkWriter
.
write
(
"poly_binary_float64.vtu"
,
Vtk
::
BINARY
,
Vtk
::
FLOAT64
);
vtkWriter
.
write
(
"poly_compressed_float64.vtu"
,
Vtk
::
COMPRESSED
,
Vtk
::
FLOAT64
);
}
Write
Preview
Supports
Markdown
0%
Try again
or
attach a new file
.
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment