[vtk-developers] how to handle FEM data in VTK?

Soeren Gebbert soerengebbert at googlemail.com
Fri Jul 4 12:35:02 EDT 2008


Hi Pierre,
i store and process FEM meshes with integration points in VTK.
I use all kind of 2d and 3d higher order elements VTK provides.

To store elements, integration points and the information about the
integration order and method,
i need to use a work-around.
The approach i have chosen and implemented stores the geometry, the integration
points and the integration Information in three different unstructured
grids. So no new data types
or cell hacks are needed. The ugrid's can be processed which the
existing vtk filter.
The three grids:

* Geometry:
An unstructured grid to store all the geometric information's. The
system nodes and elements as well as element data and point data.

* Integration points:
An unstructured grid storing the integration points as poly vertex elements.
So the integration points of one FEM element are stored as points of a
poly vertex cell
(ie: 4x5 == 20 integration points are stored as one poly vertex cell).
The order of the poly vertex cells is equal to the order of the
geometry elements.
Tensors, vectors and scalars are stored in the integration points as
point data.
You can store multiple tensors, vectors and scalars
but you can set only one of them active at a time.
Each integration point has a scalar value which defines the integration method.
Information about the method is stored in the third ugrid.
The scalar value points to the element number in the third ugrid.

* Integration order and method
The 3. unstructured grid helps to
store the information about the order and method. For each integration
order/method i
create one poly vertex cell which represent a parametric fem cell.
The integration points are stored in local parametric coordinates.
The element scalars define the number of layers and integration
points per layer and so on.

example: a model with 4 different integration types

# vtk DataFile Version 3.0
vtk output
ASCII
DATASET UNSTRUCTURED_GRID
POINTS 59 float
0.166667 0.166667 0.0469101 0.666667 0.166667 0.0469101 0.166667
0.666667 0.0469101
0.166667 0.166667 0.230765 0.666667 0.166667 0.230765 0.166667 0.666667 0.230765
0.166667 0.166667 0.5 0.666667 0.166667 0.5 0.166667 0.666667 0.5
0.166667 0.166667 0.769235 0.666667 0.166667 0.769235 0.166667 0.666667 0.769235
0.166667 0.166667 0.95309 0.666667 0.166667 0.95309 0.166667 0.666667 0.95309
0.211325 0.211325 0.0469101 0.788675 0.211325 0.0469101 0.211325
0.788675 0.0469101
0.788675 0.788675 0.0469101 0.211325 0.211325 0.230765 0.788675
0.211325 0.230765
0.211325 0.788675 0.230765 0.788675 0.788675 0.230765 0.211325 0.211325 0.5
0.788675 0.211325 0.5 0.211325 0.788675 0.5 0.788675 0.788675 0.5
0.211325 0.211325 0.769235 0.788675 0.211325 0.769235 0.211325 0.788675 0.769235
0.788675 0.788675 0.769235 0.211325 0.211325 0.95309 0.788675 0.211325 0.95309
0.211325 0.788675 0.95309 0.788675 0.788675 0.95309 0.166667 0.166667 0.0469101
0.666667 0.166667 0.0469101 0.166667 0.666667 0.0469101 0.166667
0.166667 0.230765
0.666667 0.166667 0.230765 0.166667 0.666667 0.230765 0.166667 0.166667 0.5
0.666667 0.166667 0.5 0.166667 0.666667 0.5 0.166667 0.166667 0.769235
0.666667 0.166667 0.769235 0.166667 0.666667 0.769235 0.211325
0.211325 0.0469101
0.788675 0.211325 0.0469101 0.211325 0.788675 0.0469101 0.788675
0.788675 0.0469101
0.211325 0.211325 0.230765 0.788675 0.211325 0.230765 0.211325 0.788675 0.230765
0.788675 0.788675 0.230765 0.211325 0.211325 0.5 0.788675 0.211325 0.5
0.211325 0.788675 0.5 0.788675 0.788675 0.5
CELLS 4 63
15 0 1 2 3 4 5 6 7 8 9 10 11 12 13 14
20 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34
12 35 36 37 38 39 40 41 42 43 44 45 46
12 47 48 49 50 51 52 53 54 55 56 57 58

CELL_TYPES 4
2
2
2
2

CELL_DATA 4
SCALARS NumberOfIntPoints int
LOOKUP_TABLE default
15 20 12 12
FIELD FieldData 4
Modelltype 1 4 int
23 24 15 16
NumberOfLayers 1 4 int
5 5 4 3
IntPointsLayer 1 4 int
3 4 3 4
NaturalCoordinatesDimension 1 4 int
3 3 2 2


I have implemented plenty of filters to process elements and integration
points together by using all three unstructured grids
(ie: interpolation of integration points into the system nodes or elements,
positioning of integration points and transformation of vectors and
tensors defined in the integration
points if the fem mesh is transformed and so on ...).

This is not the native way to handle FEM data, but it works very well for me.
AFAIKT this might not be the best solution. I needed to implement
many filters to handle FEM data consistently.

Best regards
Soeren

2008/7/3 Pierre JUILLARD <pierre.juillard at gmail.com>:
> Hi all,
>
> I would like to have information about the best way to handle FEM data in
> VTK file.
>
> Mainly, we use shell meshes (triangles or quadrangles) and the following
> data would be manipulated in VTK:
> - reduced element or fully integrated (1 or 4 integration point in the
> neutral fibre)
> - different integration point over thickness per element (1, 3, 5, 7, 9, 11,
> ...)
> - different integration algorithms over thickness (GAUSS, SIMPSON, ....)
> defining different weight and locations of these integration points over
> thickness
> - stress tensor at each integration point over thickness
> - strain tensor at each integration point over thickness
> - effective plastic strain at each integration point over thickness
> - thickness at each integration point in the neutral fibre
> and possibly additional results per integration points...
>
> This topic seems to have already been tackled.
> Indeed, I have read several posts about the vtkLSDynaReader class, as well
> as its documentation on
> http://public.kitware.com/VTK/doc/nightly/html/classvtkLSDynaReader.html
>
> As it is indicated in one of the post, the ultimate goal would be to benefit
> from vtk format/library to manipulate FEM data, and  from ParaView as
> PostProcessor
>
> However, from these posts, it seems that VTK is not really fitted to store
> FEM data, and for such a purpose, a good knowledge in VTK is necessary. It
> may also be necessary to extend it.
>
> Most notably, it seems currently possible to have only one scalar field,
> vector field, and tensor field per element. In the last VTK versions, has
> this "limitation" been overriden?
> Is it possible to store several fields of different types per element?
>
> Would you have some advices about the best way to proceed?
> Is there somewhere some guidelines/tutorials?
>
> It the VTK format has still some limitations, is it foreseen to override
> them in the next release?
>
> I thank you for your help.
>
> Cheers,
> _______________________________________________
> vtk-developers mailing list
> vtk-developers at vtk.org
> http://www.vtk.org/mailman/listinfo/vtk-developers
>
>



More information about the vtk-developers mailing list