Home | Show all

Skeleton data



Skeleton files

This file type is designed to store the critical points and arcs of the Morse-Smale complex (i.e. one dimensionnal filamentary structures). In skeleton files, the filamentary network is described as lists of nodes representing critical points or bifurcation points if available (see option -breakdown of skelconv), lists of segments tracing the local geometry of arcs and filaments originating from and leading to nodes and described as lists of segments. The base skeleton format is NDskl which is used internally, but this format may be converted to several other more or less complex formats of skeleton files adapted to different applications (see option -to in program skelconv, a list of available formats is displayed when running the program without argument).
nb: Within skeleton files, the segments are always oriented in the ascending direction of the arcs, which does *NOT* necessarily mean that the value of the field is necessarily increasing locally ! Indeed, the value is locally allowed to fluctuate within the persistence threshold, so the field value does not have to be strictly increasing locally along an arc.

Available formats:

  • -NDskl (Read / Write):
    This is the format of the skeleton files created with mse. It is a relatively complex binary format that contains all the information on the geometry and connectivity of the arcs of the Morse-Smale complex (i.e. filaments).


  • -NDskl_ascii (Write only):
    This ASCII format contains the same amount of information as NDskl files, but organized in a different way. In particular, filaments are not described as lists of segments but rather each filament is described by an origin node, a destination node, and a set of sampling points.


  • -segs_ascii and crits_ascii (Write only):
    This is a simplified ASCII file format designed to be easily readable but that only contains a subset of the information available in NDskl files. In particular, segs_ascii files contain a list of individual segments describing the local orientation of the arcs as well as some limited information on the filament they belong to, while crits_ascii contain information on the critical points and bifurcation points only.


  • -vtk, vtk_ascii, vtp and vtp_ascii (Write only):
    These formats are binary and ASCII legacy and XML VTK formats that are readable by several 3D visualization tools, such as VisIt or ParaView for instance.




Additional data: In addition to the topology and geometry of filamentary structures, more complete formats (i.e. NDskl, NDskl_ascii and all vtk formats) may store arbitrary additional information associated to filaments and nodes that can be used by skelconv (run skelconv filename.NDskl -info for a list of additional data available in skeleton file filename.NDskl). By default, mse stores the following additional data in skeleton type files:

  • -persistence / persistence_ratio / persistence_nsigmas (Nodes only) :
    The persistence (expressed as a difference, ratio or in number of sigmas) of the persistence pair containing the corresponding critical point. A negative or null value indicates that the node is not a critical point or that it does not belong to a persistence pair.


  • -persistence_pair (Nodes only):
    The index of the node that corresponds to the other critical point in the persistence pair (indices start at 0). The value is the index of the current node itself when it is not a critical point or it does not belong to a persistence pair.


  • -parent_index / parent_log_index (Nodes only):
    For each extrema only (i.e. minima and maxima), the index of the node that corresponds to the other extremum into which it would be merged if its persistence pair was canceled (indices start at 0). This can be used to reconstruct the tree of the hierarchy of maxima and minima. The value is -1 for non extrema critical points. The difference between the two versions is that the second (parent_log_index) is the hierarchy computed from the logarithm of the field. This second version is useful only for discrete point samples whose MS-complex is obtained from the delaunay tessellation computed with delaunay_nD. Practically, parent_log_index can be used whenever persistence pairs are cancelled in order of increasing ratio (option -nsig in mse), and parent_index whenever persistence pairs are cancelled in order of increasing difference (option -cut in mse).


  • -field_value / log_field_value (Nodes and Segments) :
    The value of the field and it logarithm. For segments, the suffix _p1 and _p2 is added to indicate which extremity of the segment it corresponds to.


  • -cell (Nodes and Segments) :
    The type and index of the cell corresponding to the node / segment in the initial network (i.e. from which the skeleton was computed). For segments, the suffix _p1 and _p2 is added to indicate which extremity of the segment it corresponds to. The value is a double precision floating number whose integer part is the index of the cell and decimal part its type. For instance, the 156th vertex (i.e. 0-cell) in the cell complex is represented as 156.0, while the 123th tetrahedron is 123.3. Note that the index of the 0-cell correpond to the index of the pixel / vertices in the original network from which the skeleton was computed.


  • -robustness / robustness_ratio (Nodes and Segments) :
    The robustness of the node / segment. Robustness is a local measure of how contrasted the critical point / filament is with respect to its local background, and it is therefore a good way to select only highly contrasted subsets of the filamentary structures (see option -trimBelow in skelconv). Note that robustness, like persistence, is defined as a difference or ratio between the value of the field in two points, so the robustness threshold has the same order of magnitude as the persistence threshold.


  • - type (Segments only) :
    The type of arc the segment belongs to. The value corresponds to the lowest critical index of the two critical points at the extremities of the arc the segments belongs to. For instance, in dimension D, the ridges (filaments) have type D-1.


NDskl format

This is the native binary format of DisPerSE. Functions for reading and writing NDskl format in C can be found within the file ${DISPERSE_SRC}/src/C/NDskeleton.c (see functions Load_NDskel, Save_NDskel and also Create_NDskel). In this format, the arcs of the Morse-Smale complex are described as arrays of nodes and segments with data associated to each of them. Each node contains information on critical points (which may also be bifurcations points or trimmed extremities of filaments). Each segment contains information on the local geometry of the arcs and global properties of the skeleton.

When using the C functions from Disperse, data is loaded into the following C structure which is close to the actual structure of the file (see file ${DISPERSE_SRC}/src/C/NDskeleton.h):

struct NDskl_seg_str {
   int nodes[2];  // index of the nodes at the extremity of the arc. Segment is oriented from nodes[0] toward nodes[1]
   float *pos;  // points to appropriate location in segpos
   int flags;  // non null if on boundary 
   int index;  // index of the segment in the Seg array
   double *data;  // points to the nsegdata supplementary data for this segment
   struct NDskl_seg_str *Next; // points to next segment in the arc, NULL if extremity (->connects to nodes[1])
   struct NDskl_seg_str *Prev; // points to previous segment in the arc, NULL if extremity (->connects to nodes[0])
 }; 
 typedef struct NDskl_seg_str NDskl_seg;
 struct NDskl_node_str {
   float *pos;  // points to appropriate location in nodepos
   int flags;  // non null if on boundary 
   int nnext;  // number of arcs connected
   int type;  // critical index
   int index;  // index of the node in the Node array
   int *nsegs; // number pf segments in each of the the nnext arcs
   double *data;   // points to the nnodedata supplementary data for this segment
   struct NDskl_node_str **Next;  // points to the node at the other extremity of each of the nnext arcs
   struct NDskl_seg_str **Seg;   // points to the first segment in the arcs, starting from current node.
 }; 
 typedef struct NDskl_node_str NDskl_node;
 typedef struct NDskel_str {
   char comment[80];
   
   int ndims;  // number of dimensions
   int *dims;  // dimension of the underlying grid, only meaningfull when computed from a regular grid
   double *x0;  // origin of the bounding box
   double *delta;  // dimension of the bounding box
    
   int nsegs;  // number of segments
   int nnodes;  // number of nodes
   
   int nsegdata;  // number of additional data for segments
   int nnodedata;  // number of additional data for nodes
   char **segdata_info;  // name of additional fields for segments
   char **nodedata_info;  // name of additional fields for nodes
   
   float *segpos;  // positions of the extremities of segments (2xndims coords for each segment)
   float *nodepos; // positions of the nodes (ndims coords for each segment)
   double *segdata; // additional data for segments (nsegs times nsegdata consecutive values)
   double *nodedata; // additional data for nodes (nnodes times nnodedata consecutive values)
   
   NDskl_seg *Seg;  // segment array (contains all segs)
   NDskl_node *Node;  // nodes array (contains all nodes)
 } NDskel;



The NDskl binary format is organized as follows (blocks are delimited by dummy variables indicating the size of the blocks for FORTRAN compatibility, but they are ignored in C):

NDskl binary format
fieldtypesizecomment
dummyint(4B) 1 for FORTRAN compatibility
tagchar(1B)16identifies the file type. Value : "NDSKEL"
dummyint(4B) 1
dummyint(4B) 1
commentchar(1B)80a comment on the file (string)
ndimsint(4B)1number of dimensions
dimsint(4B)20dimension of underlying grid (0=none, ndims values are read)
x0double(8B)20origin of bounding box (ndims first values are meaningful)
deltadouble(8B)20size of bounding box (ndims first values are meaningful)
nsegsint(4B)1number of segments
nnodesint(4B)1number of nodes
nsegdataint(4B)1number of additional data associated to segments.
nnodedataint(4B)1number of additional data associated to nodes.
dummyint(4B) 1
dummyint(4B) 1 this block is ommited if nsegdata=0
seg_data_infochar(1B)20×nsegdataname of the segment data field
dummyint(4B) 1 this block is ommited if nsegdata=0
dummyint(4B) 1 this block is ommited if nnodedata=0
node_data_infochar(1B)20×nnodedataname of the node data field
dummyint(4B) 1 this block is ommited if nnodedata=0
dummyint(4B) 1
segposfloat(4B)2×nsegs×ndimscoordinates of segment extremities
dummyint(4B) 1
dummyint(4B) 1
nodeposfloat(4B)nnodes×ndimscoordinates of nodes
dummyint(4B) 1
dummyint(4B) 1
segdatadouble(8B)nsegdata×nsegsvalue of the segments data fields
dummyint(4B) 1
dummyint(4B) 1
nodedatadouble(8B)nnodedata×nnodesvalue of the nodes data fields
dummyint(4B) 1
dummyint(4B) 1 following block is repeated for each node (×nnodes).
pos_indexint(4B)1index in the nodepos/nodedata arrays
flagsint(4B)1flags (identify boundary nodes)
nnextint(4B)1number of connected arcs
typeint(4B)1critical index (ndims+1 for bifurcations / trimmed arc axtremity)
indexint(4B)1index of this node
nsegsint(4B)nnextnumber of segments in each connected arc
* **blue section repeated for each connected arc (×nnext)
nextNodeint(4B)1index of the other node of the arc
nextSegint(4B)1index of the first segment on the arc
dummyint(4B) 1
dummyint(4B) 1 following block is repeated for each segment (×nsegs).
pos_indexint(4B)1index in the segpos/segdata arrays
nodesint(4B)2index of the 2 nodes at the endpoints of the current arc.
flagsint(4B)1flags (identify boundary nodes)
indexint(4B)1index of this segment
next_segint(4B)1index of the next segment in the node (-1 if none)
prev_segint(4B)1index of the previous segment in the node (-1 if none)
dummyint(4B) 1


NDskl_ascii format

An ASCII format that contains the same amount of information as NDskl files, but organized in a different way. In particular, filaments are not described as lists of segments but rather each filament is described by an origin node, a destination node, and a set of sampling points.

NDskl_ascii format
ANDSKEL header
ndims the number of dimensions
#comments go here OPTIONAL: should start with '#' if present (the 80 first characters are read and stored).
BBOX [x0_1 .. x0_d] [delta_1 .. delta_d] OPTIONAL: the bounding box, defined by the 'ndims' coordinates of the origin 'x0' and extent 'delta'.
[CRITICAL POINTS] Marks the beginning of the critical points section
ncrit The number of critical points (CP)
type pos[0] ... pos[ndims-1] value pairID boundary Info on the first CP: critical index, position, value, index of CP in the persistence pair, 0 if not on the boundary
nfil The number of filaments connected to this CP
destId[0] filId[0] Info on the first filament: index of the CP at the other extremity of the filament, and index of the filament (see filaments table below)
... One line for each filament connecting on the CP
destId[nfil-1] filId[nfil-1] Information on the last filament
.....
..... one blue block for each CP.
.....
[FILAMENTS] Marks the beginning of the filaments section
nfil Total number of filaments
CP1 CP2 nSamp index of the CP at the extremity of the first filament and number of sampling points
P[0][0] ... P[0][ndims-1]position of the first sampling point of first filament.
...One line for each sampling point of first filament.
P[nSamp-1][0] ... P[nSamp-1][ndims-1]Position of the last sampling point
.....
..... one blue block for each filament.
.....
[CRITICAL POINTS DATA] Marks the beginning of the CP data section
NF Number of fields associated to each CP.
CP_DATA_FIELD_NAME_1 Name of the first field
...
CP_DATA_FIELD_NAME_NF Name of the last field
val_1[0] ... val_NF[0] Value of each field for first CP
...
val_1[N_CP-1] ... val_NF[N_CP-1] Value of each field for last CP
[FILAMENTS DATA] Marks the beginning of the filaments data section
NF Number of fields associated to each sampling point of each filament.
FIL_DATA_FIELD_NAME_1 Name of the first field
...
FIL_DATA_FIELD_NAME_NF Name of the last field
val_1[0][0] ... val_NF[0][0] field values for first sampling point, first filament
...
val_1[0][nSamp[0]] ... val_NF[0][nSamp[0]] field values for last sampling point, first filament
.....
..... one blue block for each filament.
.....


segs_ascii format

This ASCII format is a very simple one. It consists in a list of individual segments describing the local orientation of the arcs as well as some limited information on the filament they belong. This may be used for doing statistical analysis of local properties of filaments such as the local orientation of a magnetic field for instance. Note that all the segments are oriented from the lower critical index extremity of the filament they belong toward the other extremity (which does NOT mean that the value is strictly increasing, as local fluctuations of amplitude the persistence threshold are still allowed ... ).
nb: you probably want to use option -breakdown of skelconv before using this format.

segs_ascii format
#arc segments Header
#U0 U1 U2 V0 V1 V2 value_U value_V type boundary Header identifying each column
#NDIMS NSEGS Number of dimensions and number of segments
U0 U1 U2 V0 V1 V2 value_U value_V type boundary The value of the fields for the first segment
.... One such line for each NSEGS segment



Notations:

  • -U_0,..., U_ndims : The coordinates of the first extremity of the segment. It is the one closest to the lowest critical index CP along the filament it belongs to.
  • -V_0,..., V_ndims : The coordinates of the second extremity of the segment. It is the one closest to the highest critical index CP along the filament it belongs to.
  • -value_U, value_V : Value of the field in U and V.
  • -type : the type of arc the segment belongs to. The type of an arc is the minimum critical index of the CP at the extremity of the filament the segment belong to.
  • -boundary : value can be 0, 1 or 2 when the segment is inside the domain, on the boundary, or outside (e.g. at infinity)


crits_ascii format

This ASCII format is a very simple one. It consists in a list of all the critical points with some useful information for each of them.

crits_ascii format
#critical points Header
#X0 X1 X2 value type pair_id boundary Header identifying each column
#NDIMS NPOINTS Number of dimensions and number of points
X0 X1 X2 value type pair_id boundary The value of the fields for the first point
.... One such line for each NPOINTS points



Notations:

  • -X_0,..., X_ndims : the coordinates of the pooint
  • -value : Value of the field.
  • -type : the critical index of the point. A value of NDIMS+1 indicates a bifurcation point or a trimmed extremity.
  • -pair_id : index of the other critical point in the persistence pair. C style convention is used, indices starting at 0, and points without perisstence pairs have their own index for pair_id.
  • -boundary : value can be 0, 1 or 2 when the point is inside the domain, on the boundary, or outside (e.g. at infinity)


vtk skeleton formats

VTK formats are developed for the Visualization Tool Kit library (VTK) and can be used for 3D visualization with software such as VisIt or ParaView. Skeletons are stored as VTK polygon data and can be output in fours different VTK formats:

  • -vtk: the legacy format
  • -vtk_ascii: ASCII version of the vtk format
  • -vtp: a more recently developed XML version of the vtk format,
  • -vtp_ascii: ASCII version of the vtp format



The specifications for these formats can be found in this PDF file. See also here and there for additional information.