DisPerSE - persistent structures identification Automatic identification of persistent structures in 2D or 3D. keywords: Morse complex, topology, peak, void, source, wall, filament, cosmic web, cosmology, structure identification. 2013-01-24T03:58:22+01:00 thierry sousbie urn:md5:f62b3bf39b9948523836119476648bbf Dotclear References urn:md5:29edbf335b4c6ddefb1ccb7dadacd273 3100-05-10T08:30:00+01:00 2012-06-03T17:26:03+02:00 thierry sousbie Contact & References <p>A few references: <br /> <br /> <ins>DisPerSE:</ins> (see therein for references on discrete Morse theory and persistence) <br /> <br /></p> <ul> <li>-Sousbie, T., "The persistent cosmic web and its filamentary structure - I. Theory and implementation", 2011, Monthly Notices of the Royal Astronomical Society, 414, 350</li> </ul> <p><br /></p> <ul> <li>-Sousbie, T., Pichon, C., Kawahara, H., "The persistent cosmic web and its filamentary structure - II. Illustrations", 2011, Monthly Notices of the Royal Astronomical Society, 414, 384</li> </ul> <p><br /> <br /> <ins>Books on persistence and discrete Morse-theory:</ins> <br /> <br /></p> <ul> <li>-Zomorodian A. J., "Topology for computing", 2009, Vol. 16 of Cambridge Monographs on Applied and Computational Mathematics, Cambridge University Press, Cambridge 2, 19, 32</li> </ul> <p><br /></p> <ul> <li>-Edelsbrunner, H. &amp; Harer, J.L., "Computational Topology", 2010</li> </ul> <p><br /></p> <ul> <li>-Gyulassy A., 2009, PhD thesis, Univ. California Berkeley 2, 8, 13</li> </ul> <p><br /></p> <ul> <li>-Forman, R., "A user’s guide to discrete morse theory", 2002</li> </ul> <p><br /> <br /></p> <p>There is also a swarm of scientific articles dealing with these topics that you can find on the net, just google it :)</p> Contact urn:md5:9894066ce7952daeb4012ab3a021e771 3000-05-10T07:52:00+01:00 2012-05-10T07:25:08+02:00 thierry sousbie Contact & References <p>For any question regarding DisPerSE or this website, please contact <strong style="color:blue"> sousbie<em>(AT)</em>iap<em>.</em>fr</strong>.</p> vtk field format urn:md5:815da7c1d79dc8027044ce918cca30b3 2430-01-06T06:56:00+01:00 2012-05-28T16:33:34+02:00 thierry sousbie Field data <p>VTK formats are developed for the Visualization Tool Kit library (<a href="http://localhost/dotclear/index.php?post/vtk-formats">VTK</a>) and can be used for 3D visualization with software such as <a href="https://wci.llnl.gov/codes/visit/">VisIt</a> or <a href="http://www.paraview.org/">ParaView</a>. Only regular grids field types can be converted to VTK (for particle type fields, convert to <em>NDnet</em> <a href="http://localhost/dotclear/index.php?post/networks">unstructured network</a> format first and then use <em><a href="http://localhost/dotclear/index.php?post/netconv">netconv</a> -to vtk</em>). Regular grids are stored as <em>VTK image data</em> and can be output in fours different VTK formats: <br /><br /></p> <ul> <li>-<strong>vtk</strong>: the legacy format</li> <li>-<strong>vtk_ascii</strong>: ASCII version of the vtk format</li> <li>-<strong>vti</strong>: a more recently developed XML version of the vtk format,</li> <li>-<strong>vti_ascii</strong>: ASCII version of the vtu format</li> </ul> <p><br /> <br /> The specifications for these formats can be found in this <a href="http://localhost/dotclear/index.php?post/www.vtk.org/VTK/img/file-formats.pdf" title="VTK formats PDF">PDF</a> file. See also <a href="http://www.cacr.caltech.edu/~slombey/asci/vtk/vtk_formats.simple.html">here</a> and <a href="http://mathema.tician.de/node/430" title="there">there</a> for additional information. <br /> <br /> <br /> <br /></p> SDL-image format urn:md5:cdd6677cf033684304583dc51bad63db 2430-01-05T06:56:00+01:00 2012-05-28T16:27:56+02:00 thierry sousbie Field data <p>This input format is used to read 2D regular grids encoded in a popular picture format readable by <a href="http://www.libsdl.org/projects/SDL_image/docs/SDL_image.html" title="SDL-image">SDL-image</a> library (jpg, png, gif, ...). When using this format as input to <a href="http://localhost/dotclear/index.php?post/mse">mse</a>, the field whose topology is computed is the luminosity L of the image, a weighted average of the RED, GREEN and BLUE components of the image: <br /> <br /></p> <pre>L = 0.2989*RED + 0.5870*GREEN + 0.1140*BLUE</pre> <p><br /> <br /> <br /> <br /> <br /> <br /> <br /> <br /> <br /> <br /></p> survey_ascii format urn:md5:0642003343a6eb7fa09f2c323576fd94 2430-01-04T06:56:00+01:00 2012-05-28T17:16:07+02:00 thierry sousbie Field data <p>This is a simple ASCII format whose main purpose is to easily read coordinates of discretely sampled astrophysical surveys (e.g. such as an <a href="http://www.sdss.org/" title="SDSS">SDSS</a> galaxy catalog). It should be considered when using spherical coordinates systems or if distance is measured by redshift for instance. <br /> <br /> In this format, particles properties are encoded in an ASCII array where each row corresponds to one particle and each column to one property. Each column must have a name defined in the header (the first line of the file, starting with character <em>#</em>). An survey_ascii file may look like this: <br /> <br /></p> <pre># ra dec z my_field +2.115401e+02 +6.313433e-01 +2.666800e-01 1 +2.115633e+02 +7.550188e-01 +1.259900e-01 0 +2.115687e+02 +8.108763e-01 +3.646600e-01 0 +2.117158e+02 +6.393598e-01 +1.143600e-01 1 +2.116826e+02 +6.528485e-01 +2.455700e-01 1 +2.116993e+02 +6.509297e-01 +1.199000e-01 0 +2.115738e+02 +7.772653e-01 +3.240600e-01 0 +2.116198e+02 +6.950604e-01 +1.987300e-01 0 +2.116773e+02 +7.085776e-01 +2.561900e-01 1</pre> <p><br /> <br /> The name of a column defines its role if it matches one of the following keywords: <br /> <br /></p> <table border="1"><caption> header keywords</caption><tr align="center"><td width="20%" style="background:grey">Keyword</td><td width="40%" style="background:grey">Meaning</td></tr> <tr align="center"><td style="background:grey"> px </td><td> X coordinate of the particle</td></tr> <tr align="center"><td style="background:grey"> py </td><td> Y coordinate of the particle</td></tr> <tr align="center"><td style="background:grey"> pz </td><td> Z coordinate of the particle</td></tr> <tr align="center"><td style="background:grey"> vx </td><td> X component of the velocity</td></tr> <tr align="center"><td style="background:grey"> vy </td><td> Y component of the velocity</td></tr> <tr align="center"><td style="background:grey"> vz </td><td> Z component of the velocity</td></tr> <tr align="center"><td style="background:grey"> id </td><td> an index associated to the particle</td></tr> <tr align="center"><td style="background:grey"> ra </td><td> The <a href="http://en.wikipedia.org/wiki/Right_ascension">right ascension</a> of the particle</td></tr> <tr align="center"><td style="background:grey"> dec </td><td> The <a href="http://en.wikipedia.org/wiki/Declination">declination</a> of the particle</td></tr> <tr align="center"><td style="background:grey"> z </td><td> The <a href="http://en.wikipedia.org/wiki/Redshift">redshift</a> of the particle.</td></tr> </table> <p><br /> <br /> When particles coordinates are defined with <em>ra</em>, <em>dec</em> or <em>z</em> (redshift), DisPerSE automatically transform them into cartesian coordinates using the standard LCDM model to compute distances (Omega_m=0.27, Omega_L=0.73 ). This transformation can be inverted on the output skeletons and networks using options <em>-toRaDecZ</em> and <em>-toRaDecDist</em> of <a href="http://localhost/dotclear/index.php?post/skelconv">skelconv</a> and <a href="http://localhost/dotclear/index.php?post/netconv">netconv</a>.</p> FITS format urn:md5:4a123cf8572bc22997acb2ed219cd137 2430-01-03T06:56:00+01:00 2012-05-28T16:21:06+02:00 thierry sousbie Field data <p>This field format is used to read 2D or 3D regular grids of arbitrary dimensions in the popular Flexible Image Transport System (<a href="http://en.wikipedia.org/wiki/FITS" title="FITSS format">FITS</a>) image format. A popular library for reading and writing FITS image in <em>C</em> or <em>FORTRAN</em> is <a href="http://heasarc.gsfc.nasa.gov/docs/software/fitsio/fitsio.html" title="FITS IO">CFITSIO</a>, which must be installed for DisPerSE to be able to use this kind of files. <br /> <br /> This file format is also used for storing <a href="http://healpix.jpl.nasa.gov/" title="HEALPIX">HEALPIX</a> tesselations of the sphere. DisPerSE will automatically detect if the FITS file contains a HEALPIX tesselation. <br /> <br /> <br /> <br /> <br /> <br /> <br /> <br /> <br /> <br /> <br /></p> NDfield format urn:md5:3aaeb963f82937bce4a9fd6addb5f958 2430-01-02T06:56:00+01:00 2012-05-28T15:56:08+02:00 thierry sousbie Field data <p>This is the native binary format of DisPerSE. Functions for reading and writing <em>NDfield</em> format in <em>C</em> can be found within the file <code>${DISPERSE_SRC}/src/C/NDfield.c</code> (see functions <em>Load_NDfield</em> and <em>Save_NDfield</em>). <br /> <br /> When using the C functions from Disperse, data is loaded into the following <em>C</em> structure which is close to the actual structure of the file (see file <code>${DISPERSE_SRC}/src/C/NDfield.h</code>):</p> <pre>#define ND_CHAR (1&lt;&lt;0) #define ND_UCHAR (1&lt;&lt;1) #define ND_SHORT (1&lt;&lt;2) #define ND_USHORT (1&lt;&lt;3) #define ND_INT (1&lt;&lt;4) #define ND_UINT (1&lt;&lt;5) #define ND_LONG (1&lt;&lt;6) #define ND_ULONG (1&lt;&lt;7) #define ND_FLOAT (1&lt;&lt;8) #define ND_DOUBLE (1&lt;&lt;9)</pre> <pre>typedef struct NDfield_str { char comment[80]; // a comment on the data int dims[NDFIELD_MAX_DIMS]; // dimensions of the grid, must be [ndims,nparticles] when data represents sample particles coordinates (i.e. when fdims_index!=0) int ndims; // number of dimensions of the space int n_dims; // number of meaningfull values in dims array int fdims_index; // if 0, the field is a regular grid of dimensions dims, else the file contains the dims[0] coordinates of dims[1] particles. int datatype; // type of the data (one of the ND_... defined above) double x0[NDFIELD_MAX_DIMS]; // origin of the bounding box double delta[NDFIELD_MAX_DIMS]; // extent of the bounding box char dummy[160]; // dummy data, for future extensions or for storing anything you want. void *val; // pointer to data long nval; // total number of particles (fdims_index==1) or pixels (fdims_index==0) int datasize; // size in bytes of datatype type. } NDfield;</pre> <p><br /> <br /> The <em>NDfield</em> binary format is organized as follows (blocks are delimited by <em>dummy</em> variables indicating the size of the blocks for FORTRAN compatibility, but they are ignored in C): <br /> <br /></p> <table style="margin: 1em auto 1em auto"><caption>NDnet binary format</caption><tr ><td width="10%" >field</td><td width="10%" >type</td><td width="10%">size</td><td width="65%" >comment</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"> for FORTRAN compatibility</td></tr> <tr valign="top"><td>tag</td><td>char(1B)</td><td>16</td><td>identifies the file type. Value : "NDFIELD"</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>ndims</td><td>int(4B)</td><td>1</td><td>number of dimensions of the embedding space</td></tr> <tr valign="top"><td>dims</td><td>int(4B)</td><td>20</td><td>size of the grid in pixels along each dimension, or [ndims,nparticles] if data represents particle coordinates (i.e. fdims_index=1)</td></tr> <tr valign="top"><td>fdims_index</td><td>int(4B)</td><td>1</td><td> 0 if data represents a regular grid, 1 if it represents coordinates of tracer particles</td></tr> <tr valign="top"><td>datatype</td><td>int(4B)</td><td>1</td><td>type of data stored (see below)</td></tr> <tr valign="top"><td>x0</td><td>double(8B)</td><td>20</td><td>origin of bounding box (first ndims val. are meaningfull)</td></tr> <tr valign="top"><td>delta</td><td>double(8B)</td><td>20</td><td>size of bounding box (first ndims val. are meaningfull)</td></tr> <tr valign="top"><td>dummy_ext</td><td>char(1B)</td><td>160</td><td>dummy data reserved for future extensions</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>data</td><td>size of datatype</td><td>N</td><td>data itself (N may be the number of pixels or ndism times the number of particles)</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> </table> <p><br /> <br /></p> <table border="1"><caption> Possible data types</caption><tr align="center"><td width="20%" style="background:grey">name</td><td width="20%" style="background:grey">size (Bytes)</td><td width="20%" style="background:grey">type</td><td width="20%" style="background:grey">value</td></tr> <tr align="center"><td style="background:grey"> ND_CHAR </td><td> 1 </td><td> integer </td><td> 1 (=1<<0)</td></tr> <tr align="center"><td style="background:grey"> ND_UCHAR </td><td> 1 </td><td> integer </td><td> 2 (=1<<1)</td></tr> <tr align="center"><td style="background:grey"> ND_SHORT </td><td> 2 </td><td> integer </td><td> 4 (=1<<2)</td></tr> <tr align="center"><td style="background:grey"> ND_USHORT </td><td> 2 </td><td> integer </td><td> 8 (=1<<3)</td></tr> <tr align="center"><td style="background:grey"> ND_INT </td><td> 4 </td><td> integer </td><td> 16 (=1<<4)</td></tr> <tr align="center"><td style="background:grey"> ND_UINT </td><td> 4 </td><td> integer </td><td> 32 (=1<<5)</td></tr> <tr align="center"><td style="background:grey"> ND_LONG </td><td> 8 </td><td> integer </td><td> 64 (=1<<6)</td></tr> <tr align="center"><td style="background:grey"> ND_ULONG </td><td> 8 </td><td> integer </td><td> 128 (=1<<7)</td></tr> <tr align="center"><td style="background:grey"> ND_FLOAT </td><td> 4 </td><td> float </td><td> 256 (=1<<8)</td></tr> <tr align="center"><td style="background:grey"> ND_DOUBLE </td><td> 8 </td><td> float </td><td> 512 (=1<<9)</td></tr> </table> NDfield_ascii format urn:md5:13c96188942b97e6721220412233a137 2430-01-02T06:56:00+01:00 2013-01-24T04:58:22+01:00 thierry sousbie Field data <p>A simple ASCII version of the <a href="http://localhost/dotclear/index.php?post/NDfield-format">NDfield</a> format used to represent either uniform grids or particles coordinates. Data type is always interpreted as double precision floating point. <br /> <br /></p> <table border="1"><caption> NDfield_ascii format</caption><tr ><td width="50%" ></td><td width="50%" ></td></tr> <tr valign="top"><th>ANDFIELD COORDS</th><td>header (COORDS keyword is optional, if present, values are interpreted as coordinates, or else, as pixel values)</td></tr> <tr valign="top"><th>[N_1 ... N_ndims]</th><td>size of the grid along each of the ndims dimensions (in pixels) or number of dimensions and number of particles if COORDS keyword is in the header (= [ndims,nparticles])</td></tr> <tr valign="top"><th>BBOX [X0_1 ... X0_ndims] [delta_1 ... delta_ndims]</th><td>OPTIONAL: bounding box definition (origin and extent)</td></tr> <tr valign="top"><th>V_1 V_2</th><td> values (any number of values can be on the same line)</td></tr> <tr valign="top"><th>V_3</th><td> no particular formating is required</td></tr> <tr valign="top"><th>...</th><td> Repeat enough time to define a number of values equal to the number of pixels or the number of particles (if COORDS keyword is in the header)</td></tr> </table> Field files urn:md5:c641a01169a789e3018403e8945394c8 2430-01-01T06:56:00+01:00 2012-07-13T14:51:29+02:00 thierry sousbie Field data <p>This file type is designed to store scalar fields sampled over regular grids (i.e. regular images in 2D, 3D and more) or, by extension, coordinates of tracer particles sampling an underlying density field. The main field format is <em>NDfield</em>, which allows to store any type of data (integer, simple or double precision floating point, ...) over arbitrary sized regular grids or as coordinates of tracer particles. <br /> <br /> <strong><ins>Available formats</ins></strong>: <br /> <br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/NDfield-format">NDfield</a></strong> (Read / Write):<br />This is the format used internally in DisPerSE, it is efficient and generic as it can store regular grids or sample particles coordinates indifferently.</li> </ul> <p><br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/NDfield_ascii-format">NDfield_ascii</a></strong> (Read / Write):<br />A simple ASCII version of the NDfield format, it is as versatile but data is always considered to be double precision floating point.</li> </ul> <p><br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/FITS-format">FITS</a></strong> (Read only):<br />Used to read 2D or 3D regular grids of arbitrary dimensions in the popular Flexible Image Transport System (<a href="http://en.wikipedia.org/wiki/FITS" title="FITSS format">FITS</a>) image format. Also used for <a href="http://healpix.jpl.nasa.gov/" title="HEALPIX">HEALPIX</a> tessellations of the sphere.</li> </ul> <p><br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/survey_ascii-format">survey_ascii</a></strong>(Read only):<br />A simple ASCII format whose main purpose is to easily read coordinates of discretely sampled astrophysical surveys (e.g. such as an <a href="http://www.sdss.org/" title="SDSS">SDSS</a> galaxy catalog). Try this if you use spherical coordinates or if distance is measured by redshift for instance.</li> </ul> <p><br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/SDL-image">SDL-image</a></strong> (Read only):<br />This input format is used to read 2D regular grids encoded in a popular picture format readable by <a href="http://www.libsdl.org/projects/SDL_image/docs/SDL_image.html" title="SDL-image">SDL-image</a> library (jpg, png, gif, ...).</li> </ul> <p><br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/vtk-field-format">vtk</a></strong>, <strong><a href="http://localhost/dotclear/index.php?post/vtk-field-format">vtk_ascii</a></strong>, <strong><a href="http://localhost/dotclear/index.php?post/vtk-field-format">vti</a></strong> and <strong><a href="http://localhost/dotclear/index.php?post/vtk-field-format">vti_ascii</a></strong> (Write only):<br />These formats are binary and ASCII legacy and XML <a href="http://www.vtk.org/">VTK</a> formats that are readable by several 3D visualization tools, such as <a href="https://wci.llnl.gov/codes/visit/">VisIt</a> or <a href="http://www.paraview.org/">ParaView</a> for instance.</li> </ul> <p><br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/networks">NDnet</a></strong> (Write only):<br /> Using this format, field files representing tracer particles coordinates can be converted to <a href="http://localhost/dotclear/index.php?post/networks">unstructured networks</a> (use <em><a href="http://localhost/dotclear/index.php?post/netconv">netconv</a></em> to convert <em>NDnet</em> files to other network formats). Regular grids conversion is not implemented yet due to the huge size of the output network.</li> </ul> <p><br /> <br /></p> vtk network format urn:md5:fb5abf70ccbb5abba106e974c9f74000 2426-01-04T13:14:00+01:00 2012-05-27T19:26:02+02:00 thierry sousbie Network data <p>VTK formats are developed for the Visualization Tool Kit library (<a href="http://localhost/dotclear/index.php?post/vtk-formats">VTK</a>) and can be used for 3D visualization with software such as <a href="https://wci.llnl.gov/codes/visit/">VisIt</a> or <a href="http://www.paraview.org/">ParaView</a>. Networks are stored as <em>VTK unstructured network data</em> and can be output in fours different VTK formats: <br /><br /></p> <ul> <li>-<strong>vtk</strong>: the legacy format</li> <li>-<strong>vtk_ascii</strong>: ASCII version of the vtk format</li> <li>-<strong>vtu</strong>: a more recently developed XML version of the vtk format,</li> <li>-<strong>vtu_ascii</strong>: ASCII version of the vtu format</li> </ul> <p><br /> <br /> The specifications for these formats can be found in this <a href="http://localhost/dotclear/index.php?post/www.vtk.org/VTK/img/file-formats.pdf" title="VTK formats PDF">PDF</a> file. See also <a href="http://www.cacr.caltech.edu/~slombey/asci/vtk/vtk_formats.simple.html">here</a> and <a href="http://mathema.tician.de/node/430" title="there">there</a> for additional information. <br /> <br /> <br /> <br /> <br /></p> ply format urn:md5:776d9154fe52b9084d1ffcd74bdebf79 2426-01-03T13:14:00+01:00 2012-05-27T19:27:27+02:00 thierry sousbie Network data <p>The PLY format is a relatively generic file format designed to store three dimensional data from 3D scanners with the possibility of associating properties to the polygons. Information on this format can be found on <a href="http://en.wikipedia.org/wiki/PLY_(file_format)" title="PLY format">wikipedia</a> (see also the <em>External links</em> section). A very efficient C library for reading and writing PLY files in ASCII or binary format is <a href="http://w3.impa.br/~diego/software/rply/" title="RPly">RPly</a> (DisPerSE uses it for PLY files I/O, see files <em>NDnet_PLY_IO.c</em> and <em>NDnet_PLY_IO.h</em>). <br /> <br /> A typical header for a PLY file readable by <a href="http://localhost/dotclear/index.php?post/netconv">netconv</a> or <a href="http://localhost/dotclear/index.php?post/mse">mse</a> is as follows:</p> <pre>ply format ascii 1.0 element bbox 1 property list uchar double x0 property list uchar double delta element vertex 77595 property float x property float y property float z property double field_value element face 482867 property list uchar uint vertex_indices end_header</pre> <p><br /> <br /></p> <ul> <li>-<strong>format</strong> may be ASCII or little / big endian binary</li> </ul> <p><br /></p> <ul> <li>-<strong>bbox</strong> element is used to define a bounding box if available (x0 is its origin and delta its extent)</li> </ul> <p><br /></p> <ul> <li>-<strong>coordinates</strong> of the vertices are given as vertex properties labeled <em>x</em>, <em>y</em>, and <em>z</em> or <em>x0</em>, <em>x1</em>, ... (the number of this properties gives the dimension of the embedding space)</li> </ul> <p><br /></p> <ul> <li>-<strong>faces</strong> are defined by the property vertex_indices, each element corresponding to a list of vertices. A cell is always supposed to be a simplex, so the number of vertices determine the dimension spanned by the complex (a 2D complex may be embedded in a 3D space).</li> </ul> <p><br /></p> <ul> <li>-<strong>additional properties</strong> may be defined for cells and vertices. In particular, a vertex property labeled <em>field_value</em> would be used as input function in <a href="http://localhost/dotclear/index.php?post/mse">mse</a>.</li> </ul> NDnet_ascii format urn:md5:e274f99d05092ef3465a29a5b43d7631 2426-01-02T13:14:00+01:00 2012-07-18T11:17:00+02:00 thierry sousbie Network data <p>This ASCII format is a simpler version of the <a href="http://localhost/dotclear/index.php?post/NDnet-format">NDnet</a> format, designed to be fully compatible but restricted to simplicial networks. It is easy to read and write and should probably be used for reasonably sized data sets. <br /> <br /> <ins>Note</ins>: The scalar function whose MS-complex is computed by <a href="http://localhost/dotclear/index.php?post/mse">mse</a> can be stored as an additional data field named 'field_value' (case sensitive).</p> <p><br /> <br /></p> <table border="1"><caption> NDnet_ascii format</caption><tr ><td width="50%" ></td><td width="50%" ></td></tr> <tr valign="top"><th> ANDNET</th><td> header</td></tr> <tr valign="top"><th> ndims</th><td> the number of dimensions</td></tr> <tr valign="top"><th> #comments go here</th><td> OPTIONAL: should start with '#' if present (the 80 first characters are read and stored).</td></tr> <tr valign="top"><th> BBOX [x0_1 .. x0_d] [delta_1 .. delta_d]</th><td> OPTIONAL: the bounding box, defined by the 'ndims' coordinates of the origin 'x0' and extent 'delta'.</td></tr> <tr valign="top"><th> nv</th><td> number of vertices</td></tr> <tr valign="top"><th> vx[0] vy[0] ...</th><td> the ndims coordinates of the first vertex</td></tr> <tr valign="top"><th> ...</th><td> One line for each vertex</td></tr> <tr valign="top"><td style="background:red" colspan="2" >Simplices definition (0-cells up to ndims-cells). One blue block should be added for each type of explicitly defined simplex. Note that only the highest dimension cells are sufficient to define a complex.</td></tr> <tr valign="top"><th style="background:lightblue"> T N</th><td style="background:lightblue"> network has N T-simplices (each T-simplex has (T+1) vertices).</td></tr> <tr valign="top"><th style="background:lightblue"> i[0] j[0] ...</th><td style="background:lightblue"> the T+1 indices (start at 0) of the vertices of the first T-simplex.</td></tr> <tr valign="top"><th style="background:lightblue"> ...</th><td style="background:lightblue"> one line for each of the N T-simplices</td></tr> <tr valign="top"><th> [ADDITIONAL_DATA]</th><td> OPTIONAL: indicate the beginning of the additional data section</td></tr> <tr valign="top"><th> additional_data_name_1</th><td> name of the additional data (e.g. field_value for mse input files)</td></tr> <tr valign="top"><th> T</th><td> type of simplex it is associated to (T-simplex, 0 means vertices)</td></tr> <tr valign="top"><th> val[0]</th><td> value for the first T-simplex</td></tr> <tr valign="top"><th> ...</th><td> One line for each T-simplex</td></tr> </table> NDnet format urn:md5:5ef6902f2833bcc58788d51806dd11ba 2426-01-01T13:14:00+01:00 2013-01-24T04:57:20+01:00 thierry sousbie Network data <p>This is the native binary format of DisPerSE. Functions for reading and writing <em>NDnet</em> format in <em>C</em> can be found within the file <code>${DISPERSE_SRC}/src/C/NDnetwork.c</code> (see functions <em>Load_NDnetwork</em> and <em>Save_NDnetwork</em>). The format may seem relatively complex, but most of it is actually optional and not used in disperse (only simplicial complexes are used in DisPerSE). To create DisPerSE input files, it is only necessary to define the highest dimensional n-simplices as a list of (n+1) vertices (see also function <em>CreateNetwork</em>). <br /> <br /> <ins>Note</ins>: The scalar function whose MS-complex is computed by <a href="http://localhost/dotclear/index.php?post/mse">mse</a> can be stored as an additional data field named 'field_value' (case sensitive). <br /> <ins>Warning</ins>: in the following, for legacy reasons, the terms <em>n-face</em> and <em>n-cell</em> are used indifferently to designate polygons of dimension <em>n</em> (which are always simplexes in DisPerSE). <br /> <br /> When using the C functions from Disperse, data is loaded into the following <em>C</em> structure which is close to the actual structure of the file (see file <code>${DISPERSE_SRC}/src/C/NDnetwork.h</code>):</p> <pre>typedef struct { int type; // the cell-type char name[255]; // name of the field double *data; // value for each of the nfaces[n] n-cells } NDnetwork_Data;</pre> <pre>// NDnetwork_SupData is not used in disperse ... typedef struct { int type; char name[255]; int datasize; char datatype[255];// a string to identity how data should be casted void *data; } NDnetwork_SupData;</pre> <pre>typedef struct { char comment[80]; int periodicity; int ndims; // the number of spatial dimensions int ndims_net; // number of dimension of the network itself (e.g. 2 for a sphere embedded in 3D) int isSimpComplex; // 1 if network is a simplicial complex (always true in disperse) double *x0; // origin of the bounding box double *delta; // size of the bounding box int indexSize; // size of NDNET_UINT type in Bytes int cumIndexSize; // size of NDNET_IDCUMT type in Bytes char dummy[160-4*2]; // dummy data reserved for future extensions NDNET_UINT nvertex; // total number of vertices float *v_coord; //vertices coodinates (X_0,Y_0,Z_0,X_1,Y_1,...,Z_nvertex-1) NDNET_UINT *nfaces; // number of cells of a given type t is given by nfaces[t] int *haveVertexFromFace; // haveVertexFromFace[n] is 1 if we have an explicit definition of the n-cells (at least one type of cell must be defined). NDNET_IDCUMT **f_numVertexIndexCum;// cumulative number of vertice in the t-cells, NULL when cells are simplexes (isSimpComplex=1) NDNET_UINT **f_vertexIndex; // list of vertices defining the n-cells is stored in f_vertexIndex[n], all vertices being enumerated for each cell (the indices of the vertices in the kth n-cell start at f_vertexIndex[n][(n+1)*k] ) // see also macro NUM_VERTEX_IN_FACE(net,type,face) and VERTEX_IN_FACE(net,type,face) //This may be computed internally within DisPerSE but does not need to be defined explicitely int *haveFaceFromVertex; // haveFaceFromVertex[n] is 1 if we have an explicit list of all the n-cells that contain each vertex (used to navigate within the network) NDNET_IDCUMT **v_numFaceIndexCum; // cumulative number of t-cells a vertex v belongs to NDNET_UINT **v_faceIndex; // indices of the t-cells in the co-boundary of v ( the list of n-cells of vertex k starts at v_faceIndex[n][net-&gt;v_numFaceIndexCum[n][k]] and ends at v_faceIndex[n][net-&gt;v_numFaceIndexCum[n][k+1]] ) // see also macro NUM_FACE_IN_VERTEX(net,type,vertex) and FACE_IN_VERTEX(net,type,vertex) // This can become extremely memory heavy ... NOT used in DisPerSE int **haveFaceFromFace; // haveFaceFromFace[k][n] is 1 if we have an explicit list of all the n-cells that have a boundary/co-boundary relation with each k-cell (used to navigate within the network) NDNET_IDCUMT ***f_numFaceIndexCum; // cumulative number of n-cells having a boundary / co-boundary relation with each k-cell: f_numFaceIndexCum[k][n] NDNET_UINT ***f_faceIndex; // indices of the cells (similar to v_faceIndex) // see also macro NUM_FACE_IN_FACE(net,ref_type,ref_face,type) and FACE_IN_FACE(net,ref_type,ref_face,type) int haveVFlags; // do we have flags associated to each vertex ? int *haveFFlags; // do we have flags associated to each n-cell ? unsigned char *v_flag; // nvertex flag values (1 for each vertex) or NULL unsigned char **f_flag; // nfaces[n] flag values (1 of each n-cell) or NULL int ndata; // number of additional data fields. NDnetwork_Data *data; // array of all additionnal data (data in total) int nsupData; NDnetwork_SupData *supData; } NDnetwork;</pre> <p><br /> <br /> The <em>NDnet</em> binary format is organized as follows (blocks are delimited by <em>dummy</em> variables indicating the size of the blocks for FORTRAN compatibility, but they are ignored in C): <br /> <br /></p> <table style="margin: 1em auto 1em auto"><caption>NDnet binary format</caption><tr valign="top"><td width="10%">field</td><td width="10%">type</td><td width="10%">size</td><td width="65%">comment</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"> for FORTRAN compatibility</td></tr> <tr valign="top"><td>tag</td><td>char(1B)</td><td>16</td><td>identifies the file type. Value : "NDNETWORK"</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>ndims</td><td>int(4B)</td><td>1</td><td>number of dimensions of the embedding space</td></tr> <tr valign="top"><td>ndims_net</td><td>int(4B)</td><td>1</td><td>ndims spanned by the network (=ndims by default)</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>comment</td><td>char(1B)</td><td>80</td><td>a comment on the file (string)</td></tr> <tr valign="top"><td> periodicity</td><td>int(4B)</td><td>1</td><td> 0=non periodic, if p^th bit is set, boundary are periodic along dimension p</td></tr> <tr valign="top"><td>isSimpComplex</td><td>int(4B)</td><td>1</td><td> 1 if network is made of simplices (must be 1 for DisPerSE)</td></tr> <tr valign="top"><td>x0</td><td>double(8B)</td><td>ndims</td><td>origin of bounding box</td></tr> <tr valign="top"><td>delta</td><td>double(8B)</td><td>ndims</td><td>size of bounding box</td></tr> <tr valign="top"><td>index_size</td><td>int(4B)</td><td>1</td><td>size of NDNET_UINT integer format in Bytes</td></tr> <tr valign="top"><td>cumindex_size</td><td>int(4B)</td><td>1</td><td>size of NDNET_IDCUMT integer format in Bytes</td></tr> <tr valign="top"><td>dummy_ext</td><td>char(1B)</td><td>152</td><td>dummy data reserved for future extensions</td></tr> <tr valign="top"><td>nvertex</td><td>NDNET_UINT</td><td>1</td><td>number of vertices</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>v_coords</td><td>float(4B)</td><td>ndims&timesnvertex</td><td>coordinates of the vertices [X0,Y0, ...]</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>nfaces</td><td>NDNET_UINT</td><td>ndims+1</td><td>number of cells of each type (N0,N1,...)</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>haveVertexFromFace</td><td>int(4B)</td><td>ndims+1</td><td>are n-cells explicitly defined ? (0=no, 1=yes)</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">*</td><td style="background:grey">*</td><td style="background:grey">*</td><td style="background:red"> next 3 lines are repeated for each (ndims+1) possible cells type , only if haveVertexFromFace[n] is true.</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> <tr valign="top"><td>f_vertexIndex[n]</td><td>NDNET_UINT</td><td>(n+1)&timesnfaces[n]</td><td style="background:lightblue">list of (n+1) vertex indices for each n-cell</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>haveFaceFromVertex</td><td>int(4B)</td><td>ndims+1</td><td>are n-cells in the co-boundary of each vertex explicitly defined ? (0=no, 1=yes)</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">*</td><td style="background:grey">*</td><td style="background:grey">*</td><td style="background:red"> next 6 lines are repeated for each (ndims+1) possible cells type, only if haveFaceFromVertex[n] is true.</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> <tr valign="top"><td>numFaceIndexCum[n]</td><td>NDNET_IDCUMT</td><td>nvertex+1</td><td style="background:lightblue">cumulative count of n-cells on vertices co-boundary</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> <tr valign="top"><td>v_faceIndex[n]</td><td>NDNET_UINT</td><td>numFaceIndexCum[n]</td><td style="background:lightblue">list of n-cells on the co-boundary of each vertex (vertex i have numFaceIndexCum[n][i+1]-numFaceIndexCum[n][i] of them)</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>haveFaceFromFace</td><td>int(4B)</td><td>(ndims+1)^2</td><td>are n-cells in the co-boundary of each k-cell explicitly defined ? (0=no, 1=yes)</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>*</td><td>*</td><td>*</td><td style="background:red"> This section describes boundary relation between n-cells and k-cells. It is usually empty in DisPerSE (see NDnetwork.c for details) so SKIP IT :)</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>haveVFlags</td><td>int(4B)</td><td>1</td><td> 1 if vertex flags are defined</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>*</td><td>*</td><td>*</td><td style="background:red"> next 3 lines are skipped if haveVFlags=0</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> <tr valign="top"><td>v_flag</td><td>uchar(1B)</td><td>nvertex</td><td style="background:lightblue"> value of the flags for vertices</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>haveFFlags</td><td>int(4B)</td><td>ndims+1</td><td> 1 if flags are defined for n-cells</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>*</td><td>*</td><td>*</td><td style="background:red"> next 3 lines are repeated for each n-cell such that haveFFlags[n]=1</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> <tr valign="top"><td>f_flag[n]</td><td>uchar(1B)</td><td>nfaces[n]</td><td style="background:lightblue"> value of the flags for n-cells</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>ndata</td><td>int(4B)</td><td>1</td><td> total number of additional fields</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>*</td><td>*</td><td>*</td><td style="background:red"> next 7 lines are repeated for each additional field (&timesndata)</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> <tr valign="top"><td>type</td><td>int(4B)</td><td>1</td><td style="background:lightblue"> the type of cells (0=vertex, n = n-simplex)</td></tr> <tr valign="top"><td>name</td><td>char(1B)</td><td>255</td><td style="background:lightblue"> name of the supplementary data</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> <tr valign="top"><td>data</td><td>double(8B)</td><td>N</td><td style="background:lightblue"> data associated to cells or vertices. N is nfaces[n] or nvertex depending on the type value.</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:lightblue"></td></tr> </table> Network files urn:md5:155e7207695dd84931742192590c0813 2425-05-09T06:54:00+01:00 2012-07-18T11:18:40+02:00 thierry sousbie Network data <p>This file type is designed to store any kind of unstructured network (such as delaunay tesselations or more generally cell complexes). In DisPerSE, its usage is restricted to networks of simplices though, and it is mainly used to store ascending and descending manifolds (voids and walls for instance) and persistence pairs as output by <a href="http://localhost/dotclear/index.php?post/mse">mse</a> or Delaunay tessellations as output by <a href="http://localhost/dotclear/index.php?post/Usage">delaunay_nD</a> (<a href="http://localhost/dotclear/index.php?post/skeleton-formats">skeleton files</a> can also be converted to networks using <em><a href="http://localhost/dotclear/index.php?post/skelconv">skelconv</a> -to NDnet</em>). Within network files, networks are represented by setz of vertices and cells of any dimension. A n-cell is a cell of dimension <em>n</em>, which is therefore described, in the case of a simplicial network, as a set of <em>n+1</em> vertices index. In the case of a simplicial <em>complex</em>, only the highest dimensional cell have to be explicitly given, but other type of cells may also be specified. Indeed, extended manifolds for instance are not described as complexes: an ascending 0-manifold in 3D is a set of tetrahedrons (3-cells), triangles (representing ascending 1-manifolds on its boundary), segments (ascending 2-manifolds on its boundary) and vertices (ascending 3-manifolds / critical points). Note that additional information can also be associated to each type of cell (see below). <br /> The base network format is <a href="http://localhost/dotclear/index.php?post/NDnet-format">NDnet</a> which is used internally, but this format may be converted to several other more or less complex formats of network files adapted to different applications (see option <em>-to</em> in program <a href="http://localhost/dotclear/index.php?post/netconv">netconv</a>, a list of available formats is displayed when running the program without argument).<br /> <br /> <br /> <strong><ins>Available formats</ins></strong>: <br /> <br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/NDnet-format">NDnet</a></strong> (Read / Write):<br /> This is the format of the network files created or red by <a href="http://localhost/dotclear/index.php?post/mse">mse</a>. It is a relatively complex binary format (it is actually more complex than needed as it is designed to store generic non-simplicial networks) that contains all the information on the geometry and topology of unstructured networks as well as additional data associated to each type of cells.</li> </ul> <p><br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/NDnet_ascii-format">NDnet_ascii</a></strong> (Read / Write):<br /> This ASCII format contains the same amount of information as <em>NDnet</em> files, but restricted to simplicial networks. It is easy to read and write so it may be used to write reasonably sized networks used as input for <em><a href="http://localhost/dotclear/index.php?post/mse">mse</a></em>.</li> </ul> <p><br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/ply-format">PLY</a></strong> and <strong><a href="http://localhost/dotclear/index.php?post/ply-format">PLY_ascii</a></strong> (Read binary only / Write):<br /> This is a rather popular and simple binary or ASCII format that can be used as an interface with other software.</li> </ul> <p><br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/vtk-network-format">vtk</a></strong>, <strong><a href="http://localhost/dotclear/index.php?post/vtk-network-format">vtk_ascii</a></strong>, <strong><a href="http://localhost/dotclear/index.php?post/vtk-network-format">vtu</a></strong> and <strong><a href="http://localhost/dotclear/index.php?post/vtk-network-formats">vtu_ascii</a></strong> (Write only):<br /> These formats are binary and ASCII legacy and XML <a href="http://www.vtk.org/">VTK</a> formats that are readable by several 3D visualization tools, such as <a href="https://wci.llnl.gov/codes/visit/">VisIt</a> or <a href="http://www.paraview.org/">ParaView</a> for instance.</li> </ul> <p><br /> <br /> <strong><ins>Additional data</ins></strong>: In addition to the topology and geometry of the network, arbitrary additional information may be associated to each type of cell. Run <em><a href="http://localhost/dotclear/index.php?post/netconv">netconv</a> filename -info</em> for a list of additional data available in the network file <em>filename</em>. By default, the name of additional data added by <a href="http://localhost/dotclear/index.php?post/mse">mse</a> is relatively explicit, it includes (see also the <em>additional data</em> section of the <a href="http://localhost/dotclear/index.php?post/skeleton-formats">skeleton file format</a> description): <br /> <br /></p> <ul> <li>-<strong>field_value</strong> / <strong>log_field_value</strong> :<br /> The value of the field and its logarithm. The tag <em>field_value</em> corresponds to the input function for <em><a href="http://localhost/dotclear/index.php?post/mse">mse</a></em>, whose Morse-Smale complex is to be computed.</li> </ul> <p><br /></p> <ul> <li>-<strong>cell</strong>:<br /> The <em>type</em> and <em>index</em> of a cell in the original network (prefix may be added). 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.</li> </ul> <p><br /></p> <ul> <li>-<strong>type</strong>:<br /> This usually corresponds to the critical index of a critical point (for instance, vertices of persistence pairs networks), or the type of a persistence pair (i.e. the minimum critical index of the CP in the pair, for segments of persistence pairs networks).</li> </ul> <p><br /></p> <ul> <li>-<strong>index</strong><br /> Usually the index of a vertex (e.g. for persistence pairs, additional segment data tagged <em>up_index</em> and <em>down_index</em> correspond to the indices of the vertices with lowest and highest critical index in the persistence pair respectively).</li> </ul> <p><br /></p> <ul> <li>-<strong>persistence</strong> / <strong>persistence_ratio</strong> / <strong>persistence_nsigmas</strong> :<br /> The persistence (expressed as a difference, ratio or in <em>number of sigmas</em>) of the persistence pair containing the corresponding critical point. A negative or null value indicates that <em>persistence</em> is not relevant to this particular cell.</li> </ul> <p><br /><a name="parent"></a></p> <ul> <li>-<strong>parent_index</strong> / <strong>parent_log_index</strong> (vertices only):<br /> For persistence pairs type networks, for each vertex representing an extremum (i.e. minima and maxima), the index of the vertex 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 (<em>parent_log_index</em>) is the hierarchy computed from the logarithm of the field. The second version is useful only for discrete point samples whose MS-complex is obtained from the delaunay tessellation computed with <em><a href="http://localhost/dotclear/index.php?post/Usage">delaunay_nD</a></em>. Practically, <em>parent_log_index</em> can be used whenever persistence pairs are cancelled in order of increasing ratio (option <em><a href="http://localhost/dotclear/index.php?post/mse#nsig">-nsig</a></em> in <em>mse</em>), and <em>parent_index</em> whenever persistence pairs are cancelled in order of increasing difference (option <em><a href="http://localhost/dotclear/index.php?post/mse#cut">-cut</a></em> in <em>mse</em>).</li> </ul> <p><br /><a name="source"></a></p> <ul> <li>-<strong>source_cell</strong> / <strong>source_index</strong>:<br /> For networks representing manifolds (voids, walls, ... obtained with option <a href="http://localhost/dotclear/index.php?post/mse#dumpmanifolds">-dumpManifolds</a> of mse), this represents for each simplex the critical point from which the manifold it belongs to originates (for instance, the minimum corresponding to a void, or the saddle point corresponding to a filament). In <em>source_cell</em>, the critical points is represented by its cell in the initial cell complex (see <strong>cell</strong> above), while <em>source_index</em> gives the index of the critical point in the <a href="http://localhost/dotclear/index.php?post/skeleton-formats">skeleton</a> file or persistence pair network obtained with <a href="http://localhost/dotclear/index.php?post/mse">mse</a> (options <em><a href="http://localhost/dotclear/index.php?post/mse#dumparcs">-dumpArcs</a></em> and <em><a href="http://localhost/dotclear/index.php?post/mse#ppairs">-ppairs</a></em>). See also <a href="http://localhost/dotclear/index.php?post/Example-1#voids">here</a> and <a href="http://localhost/dotclear/index.php?post/Example-2#walls">there</a> in the tutorial section.</li> </ul> vtk skeleton formats urn:md5:611138afc79e8da4f1ebda37f6c51ad2 2420-05-07T15:06:00+01:00 2012-05-26T11:51:24+02:00 thierry sousbie Skeleton data <p>VTK formats are developed for the Visualization Tool Kit library (<a href="http://localhost/dotclear/index.php?post/vtk-formats">VTK</a>) and can be used for 3D visualization with software such as <a href="https://wci.llnl.gov/codes/visit/">VisIt</a> or <a href="http://www.paraview.org/">ParaView</a>. Skeletons are stored as <em>VTK polygon data</em> and can be output in fours different VTK formats: <br /><br /></p> <ul> <li>-<strong>vtk</strong>: the legacy format</li> <li>-<strong>vtk_ascii</strong>: ASCII version of the vtk format</li> <li>-<strong>vtp</strong>: a more recently developed XML version of the vtk format,</li> <li>-<strong>vtp_ascii</strong>: ASCII version of the vtp format</li> </ul> <p><br /> <br /> The specifications for these formats can be found in this <a href="http://localhost/dotclear/index.php?post/www.vtk.org/VTK/img/file-formats.pdf" title="VTK formats PDF">PDF</a> file. See also <a href="http://www.cacr.caltech.edu/~slombey/asci/vtk/vtk_formats.simple.html">here</a> and <a href="http://mathema.tician.de/node/430" title="there">there</a> for additional information. <br /> <br /> <br /> <br /> <br /></p> crits_ascii format urn:md5:a36d6123b06f743a8cfd629740502e4d 2420-05-06T15:03:00+01:00 2012-05-24T16:31:02+02:00 thierry sousbie Skeleton data <p>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. <br /> <br /></p> <table border="1"><caption> crits_ascii format</caption><tr ><td width="50%" ></td><td width="50%" ></td></tr> <tr valign="top"><th> #critical points</th><td> Header</td></tr> <tr valign="top"><th> #X0 X1 X2 value type pair_id boundary</th><td> Header identifying each column</td></tr> <tr valign="top"><th> #NDIMS NPOINTS</th><td> Number of dimensions and number of points</td></tr> <tr valign="top"><th> X0 X1 X2 value type pair_id boundary</th><td> The value of the fields for the first point</td></tr> <tr valign="top"><th> ....</th><td> One such line for each NPOINTS points</td></tr> </table> <p><br /> <br /> <ins>Notations</ins>:</p> <ul> <li>-<strong>X_0</strong>,..., <strong>X_ndims</strong> : the coordinates of the pooint</li> <li>-<strong>value</strong> : Value of the field.</li> <li>-<strong>type</strong> : the critical index of the point. A value of NDIMS+1 indicates a bifurcation point or a trimmed extremity.</li> <li>-<strong>pair_id</strong> : 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 <em>pair_id</em>.</li> <li>-<strong>boundary</strong> : value can be <em>0</em>, <em>1</em> or <em>2</em> when the point is inside the domain, on the boundary, or outside (e.g. at infinity)</li> </ul> segs_ascii format urn:md5:364e501555e6d066eda82436e0a95325 2420-05-04T15:03:00+01:00 2013-01-24T04:56:07+01:00 thierry sousbie Skeleton data <p>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 ... ). <br /> <ins><strong>nb</strong></ins>: you probably want to use option <em><a href="http://localhost/dotclear/index.php?post/skelconv#breakdown">-breakdown</a></em> of <em><a href="http://localhost/dotclear/index.php?post/skelconv">skelconv</a></em> before using this format. <br /> <br /></p> <table border="1"><caption> segs_ascii format</caption><tr ><td width="50%" ></td><td width="50%" ></td></tr> <tr valign="top"><th> #arc segments</th><td> Header</td></tr> <tr valign="top"><th> #U0 U1 U2 V0 V1 V2 value_U value_V type boundary</th><td> Header identifying each column</td></tr> <tr valign="top"><th> #NDIMS NSEGS</th><td> Number of dimensions and number of segments</td></tr> <tr valign="top"><th> U0 U1 U2 V0 V1 V2 value_U value_V type boundary</th><td> The value of the fields for the first segment</td></tr> <tr valign="top"><th> ....</th><td> One such line for each NSEGS segment</td></tr> </table> <p><br /> <br /> <ins>Notations</ins>:</p> <ul> <li>-<strong>U_0</strong>,..., <strong>U_ndims</strong> : 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.</li> <li>-<strong>V_0</strong>,..., <strong>V_ndims</strong> : 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.</li> <li>-<strong>value_U</strong>, <strong>value_V</strong> : Value of the field in U and V.</li> <li>-<strong>type</strong> : 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.</li> <li>-<strong>boundary</strong> : value can be <em>0</em>, <em>1</em> or <em>2</em> when the segment is inside the domain, on the boundary, or outside (e.g. at infinity)</li> </ul> NDskl_ascii format urn:md5:6f0d89c6789bb8662825d1d4061e3b81 2420-05-03T15:02:00+01:00 2012-05-27T15:50:00+02:00 thierry sousbie Skeleton data <p>An ASCII format that contains the same amount of information as <a href="http://localhost/dotclear/index.php?post/NDskl-format">NDskl</a> 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. <br /> <br /></p> <table border="1"><caption> NDskl_ascii format</caption><tr ><td width="50%" ></td><td width="50%" ></td></tr> <tr valign="top"><th> ANDSKEL</th><td> header</td></tr> <tr valign="top"><th> ndims</th><td> the number of dimensions</td></tr> <tr valign="top"><th> #comments go here</th><td> OPTIONAL: should start with '#' if present (the 80 first characters are read and stored).</td></tr> <tr valign="top"><th> BBOX [x0_1 .. x0_d] [delta_1 .. delta_d]</th><td> OPTIONAL: the bounding box, defined by the 'ndims' coordinates of the origin 'x0' and extent 'delta'.</td></tr> <tr valign="top"><th> [CRITICAL POINTS]</th><td> Marks the beginning of the critical points section</td></tr> <tr valign="top"><th> ncrit</th><td> The number of critical points (CP)</td></tr> <tr valign="top"><th> type pos[0] ... pos[ndims-1] value pairID boundary</th><td style="background:lightblue"> Info on the first CP: critical index, position, value, index of CP in the persistence pair, 0 if not on the boundary</td></tr> <tr valign="top"><th> nfil</th><td style="background:lightblue"> The number of filaments connected to this CP</td></tr> <tr valign="top"><th> destId[0] filId[0]</th><td style="background:lightblue"> 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)</td></tr> <tr valign="top"><th> ...</th><td style="background:lightblue"> One line for each filament connecting on the CP</td></tr> <tr valign="top"><th> destId[nfil-1] filId[nfil-1]</th><td style="background:lightblue"> Information on the last filament</td></tr> <tr valign="top"><th style="background:red">.....</th><td style="background:red"></td></tr> <tr valign="top"><th style="background:red"> .....</th><td style="background:red"> one blue block for each CP.</td></tr> <tr valign="top"><th style="background:red"> .....</th><td style="background:red"></td></tr> <tr valign="top"><th> [FILAMENTS]</th><td> Marks the beginning of the filaments section</td></tr> <tr valign="top"><th> nfil</th><td> Total number of filaments</td></tr> <tr valign="top"><th> CP1 CP2 nSamp</th><td style="background:lightblue"> index of the CP at the extremity of the first filament and number of sampling points</td></tr> <tr valign="top"><th> P[0][0] ... P[0][ndims-1]</th><td style="background:lightblue">position of the first sampling point of first filament.</td></tr> <tr valign="top"><th> ...</th><td style="background:lightblue">One line for each sampling point of first filament.</td></tr> <tr valign="top"><th> P[nSamp-1][0] ... P[nSamp-1][ndims-1]</th><td style="background:lightblue">Position of the last sampling point</td></tr> <tr valign="top"><th style="background:red">.....</th><td style="background:red"></td></tr> <tr valign="top"><th style="background:red"> .....</th><td style="background:red"> one blue block for each filament.</td></tr> <tr valign="top"><th style="background:red"> .....</th><td style="background:red"></td></tr> <tr valign="top"><th> [CRITICAL POINTS DATA]</th><td> Marks the beginning of the CP data section</td></tr> <tr valign="top"><th> NF</th><td> Number of fields associated to each CP.</td></tr> <tr valign="top"><th> CP_DATA_FIELD_NAME_1</th><td> Name of the first field</td></tr> <tr valign="top"><th> ...</th><td></td></tr> <tr valign="top"><th> CP_DATA_FIELD_NAME_NF</th><td> Name of the last field</td></tr> <tr valign="top"><th> val_1[0] ... val_NF[0]</th><td> Value of each field for first CP</td></tr> <tr valign="top"><th> ...</th><td></td></tr> <tr valign="top"><th> val_1[N_CP-1] ... val_NF[N_CP-1]</th><td> Value of each field for last CP</td></tr> <tr valign="top"><th> [FILAMENTS DATA]</th><td> Marks the beginning of the filaments data section</td></tr> <tr valign="top"><th> NF</th><td> Number of fields associated to each sampling point of each filament.</td></tr> <tr valign="top"><th> FIL_DATA_FIELD_NAME_1</th><td> Name of the first field</td></tr> <tr valign="top"><th> ...</th><td></td></tr> <tr valign="top"><th> FIL_DATA_FIELD_NAME_NF</th><td> Name of the last field</td></tr> <tr valign="top"><th> val_1[0][0] ... val_NF[0][0]</th><td style="background:lightblue"> field values for first sampling point, first filament</td></tr> <tr valign="top"><th> ...</th><td style="background:lightblue"></td></tr> <tr valign="top"><th> val_1[0][nSamp[0]] ... val_NF[0][nSamp[0]]</th><td style="background:lightblue"> field values for last sampling point, first filament</td></tr> <tr valign="top"><th style="background:red">.....</th><td style="background:red"></td></tr> <tr valign="top"><th style="background:red"> .....</th><td style="background:red"> one blue block for each filament.</td></tr> <tr valign="top"><th style="background:red"> .....</th><td style="background:red"></td></tr> </table> NDskl format urn:md5:c09ea460592168c16755715666249515 2420-05-02T15:01:00+01:00 2012-05-27T12:55:36+02:00 thierry sousbie Skeleton data <p>This is the native binary format of DisPerSE. Functions for reading and writing <em>NDskl</em> format in <em>C</em> can be found within the file <code>${DISPERSE_SRC}/src/C/NDskeleton.c</code> (see functions <em>Load_NDskel</em>, <em>Save_NDskel</em> and also <em>Create_NDskel</em>). 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. <br /> <br /> When using the C functions from Disperse, data is loaded into the following <em>C</em> structure which is close to the actual structure of the file (see file <code>${DISPERSE_SRC}/src/C/NDskeleton.h</code>):</p> <pre>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 (-&gt;connects to nodes[1]) struct NDskl_seg_str *Prev; // points to previous segment in the arc, NULL if extremity (-&gt;connects to nodes[0]) }; typedef struct NDskl_seg_str NDskl_seg;</pre> <pre> 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;</pre> <pre> 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;</pre> <p><br /> <br /> The <em>NDskl</em> binary format is organized as follows (blocks are delimited by <em>dummy</em> variables indicating the size of the blocks for FORTRAN compatibility, but they are ignored in C): <br /> <br /></p> <table style="margin: 1em auto 1em auto"><caption>NDskl binary format</caption><tr valign="top"><td width="10%">field</td><td width="10%">type</td><td width="10%">size</td><td width="65%">comment</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"> for FORTRAN compatibility</td></tr> <tr valign="top"><td>tag</td><td>char(1B)</td><td>16</td><td>identifies the file type. Value : "NDSKEL"</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>comment</td><td>char(1B)</td><td>80</td><td>a comment on the file (string)</td></tr> <tr valign="top"><td>ndims</td><td>int(4B)</td><td>1</td><td>number of dimensions</td></tr> <tr valign="top"><td>dims</td><td>int(4B)</td><td>20</td><td>dimension of underlying grid (0=none, ndims values are read)</td></tr> <tr valign="top"><td>x0</td><td>double(8B)</td><td>20</td><td>origin of bounding box (ndims first values are meaningful)</td></tr> <tr valign="top"><td>delta</td><td>double(8B)</td><td>20</td><td>size of bounding box (ndims first values are meaningful)</td></tr> <tr valign="top"><td>nsegs</td><td>int(4B)</td><td>1</td><td>number of segments</td></tr> <tr valign="top"><td>nnodes</td><td>int(4B)</td><td>1</td><td>number of nodes</td></tr> <tr valign="top"><td>nsegdata</td><td>int(4B)</td><td>1</td><td>number of additional data associated to segments.</td></tr> <tr valign="top"><td>nnodedata</td><td>int(4B)</td><td>1</td><td>number of additional data associated to nodes.</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"> this block is ommited if nsegdata=0</td></tr> <tr valign="top"><td>seg_data_info</td><td>char(1B)</td><td>20&timesnsegdata</td><td>name of the segment data field</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"> this block is ommited if nsegdata=0</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"> this block is ommited if nnodedata=0</td></tr> <tr valign="top"><td>node_data_info</td><td>char(1B)</td><td>20&timesnnodedata</td><td>name of the node data field</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"> this block is ommited if nnodedata=0</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>segpos</td><td>float(4B)</td><td>2&timesnsegs&timesndims</td><td>coordinates of segment extremities</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>nodepos</td><td>float(4B)</td><td>nnodes&timesndims</td><td>coordinates of nodes</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>segdata</td><td>double(8B)</td><td>nsegdata&timesnsegs</td><td>value of the segments data fields</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td>nodedata</td><td>double(8B)</td><td>nnodedata&timesnnodes</td><td>value of the nodes data fields</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:red"> following block is repeated for each node (&timesnnodes).</td></tr> <tr valign="top"><td>pos_index</td><td>int(4B)</td><td>1</td><td>index in the nodepos/nodedata arrays</td></tr> <tr valign="top"><td>flags</td><td>int(4B)</td><td>1</td><td>flags (identify boundary nodes)</td></tr> <tr valign="top"><td>nnext</td><td>int(4B)</td><td>1</td><td>number of connected arcs</td></tr> <tr valign="top"><td>type</td><td>int(4B)</td><td>1</td><td>critical index (ndims+1 for bifurcations / trimmed arc axtremity)</td></tr> <tr valign="top"><td>index</td><td>int(4B)</td><td>1</td><td>index of this node</td></tr> <tr valign="top"><td>nsegs</td><td>int(4B)</td><td>nnext</td><td>number of segments in each connected arc</td></tr> <tr valign="top"><td> * </td><td>*</td><td>*</td><td style="background:red">blue section repeated for each connected arc (&timesnnext)</td></tr> <tr valign="top"><td>nextNode</td><td>int(4B)</td><td>1</td><td style="background:lightblue">index of the other node of the arc</td></tr> <tr valign="top"><td>nextSeg</td><td>int(4B)</td><td>1</td><td style="background:lightblue">index of the first segment on the arc</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:red"> following block is repeated for each segment (&timesnsegs).</td></tr> <tr valign="top"><td>pos_index</td><td>int(4B)</td><td>1</td><td>index in the segpos/segdata arrays</td></tr> <tr valign="top"><td>nodes</td><td>int(4B)</td><td>2</td><td>index of the 2 nodes at the endpoints of the current arc.</td></tr> <tr valign="top"><td>flags</td><td>int(4B)</td><td>1</td><td>flags (identify boundary nodes)</td></tr> <tr valign="top"><td>index</td><td>int(4B)</td><td>1</td><td>index of this segment</td></tr> <tr valign="top"><td>next_seg</td><td>int(4B)</td><td>1</td><td>index of the next segment in the node (-1 if none)</td></tr> <tr valign="top"><td>prev_seg</td><td>int(4B)</td><td>1</td><td>index of the previous segment in the node (-1 if none)</td></tr> <tr valign="top"><td style="background:grey">dummy</td><td style="background:grey">int(4B) </td><td style="background:grey">1</td><td style="background:grey"></td></tr> <tr valign="top"></tr> </table> Skeleton files urn:md5:5a11e7d49bef99b40d28c01c537d6cad 2420-05-01T06:56:00+01:00 2012-07-18T11:10:05+02:00 thierry sousbie Skeleton data <p>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 <em><a href="http://localhost/dotclear/index.php?post/skelconv#breakdown">-breakdown</a></em> of <a href="http://localhost/dotclear/index.php?post/skelconv">skelconv</a>), 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 <a href="http://localhost/dotclear/index.php?post/NDskl-format">NDskl</a> 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 <em><a href="http://localhost/dotclear/index.php?post/skelconv#to">-to</a></em> in program <a href="http://localhost/dotclear/index.php?post/skelconv">skelconv</a>, a list of available formats is displayed when running the program without argument).<br /> <ins>nb</ins>: 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. <br /> <br /> <strong><ins>Available formats</ins></strong>: <br /> <br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/NDskl-format">NDskl</a></strong> (Read / Write):<br /> This is the format of the skeleton files created with <a href="http://localhost/dotclear/index.php?post/mse">mse</a>. 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).</li> </ul> <p><br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/NDskl_ascii-format">NDskl_ascii</a></strong> (Write only):<br /> This ASCII format contains the same amount of information as <em>NDskl</em> 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.</li> </ul> <p><br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/segs_ascii-format">segs_ascii</a></strong> and <strong><a href="http://localhost/dotclear/index.php?post/crits_ascii-format">crits_ascii</a></strong> (Write only):<br /> This is a simplified ASCII file format designed to be easily readable but that only contains a subset of the information available in <em>NDskl</em> files. In particular, <em>segs_ascii</em> 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 <em>crits_ascii</em> contain information on the critical points and bifurcation points only.</li> </ul> <p><br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/vtk-formats">vtk</a></strong>, <strong><a href="http://localhost/dotclear/index.php?post/vtk-formats">vtk_ascii</a></strong>, <strong><a href="http://localhost/dotclear/index.php?post/vtk-formats">vtp</a></strong> and <strong><a href="http://localhost/dotclear/index.php?post/vtk-formats">vtp_ascii</a></strong> (Write only):<br /> These formats are binary and ASCII legacy and XML <a href="http://www.vtk.org/">VTK</a> formats that are readable by several 3D visualization tools, such as <a href="https://wci.llnl.gov/codes/visit/">VisIt</a> or <a href="http://www.paraview.org/">ParaView</a> for instance.</li> </ul> <p><br /></p> <ul> <li>-<strong><a href="http://localhost/dotclear/index.php?post/networks">NDnet</a></strong> (Write only):<br /> Using this format, skeleton files can be converted to <a href="http://localhost/dotclear/index.php?post/networks">unstructured networks</a> (use <em><a href="http://localhost/dotclear/index.php?post/netconv">netconv</a></em> to convert <em>NDnet</em> files to other network formats).</li> </ul> <p><br /> <br /> <strong><ins>Additional data</ins></strong>: In addition to the topology and geometry of filamentary structures, more complete formats (i.e. <a href="http://localhost/dotclear/index.php?post/NDskl-format">NDskl</a>, <a href="http://localhost/dotclear/index.php?post/NDskl_ascii-format">NDskl_ascii</a> and all <a href="http://localhost/dotclear/index.php?post/vtk-formats">vtk</a> formats) may store arbitrary additional information associated to filaments and nodes that can be used by <em><a href="http://localhost/dotclear/index.php?post/skelconv">skelconv</a></em> (run <em>skelconv filename.NDskl -info</em> for a list of additional data available in skeleton file <em>filename.NDskl</em>). By default, <a href="http://localhost/dotclear/index.php?post/mse">mse</a> stores the following additional data in skeleton type files: <br /> <br /></p> <ul> <li>-<strong>persistence</strong> / <strong>persistence_ratio</strong> / <strong>persistence_nsigmas</strong> (Nodes only) :<br /> The persistence (expressed as a difference, ratio or in <em>number of sigmas</em>) 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.</li> </ul> <p><br /></p> <ul> <li>-<strong>persistence_pair</strong> (Nodes only):<br /> 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.</li> </ul> <p><br /><a name="parent"></a></p> <ul> <li>-<strong>parent_index</strong> / <strong>parent_log_index</strong> (Nodes only):<br /> 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 (<em>parent_log_index</em>) 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 <em><a href="http://localhost/dotclear/index.php?post/Usage">delaunay_nD</a></em>. Practically, <em>parent_log_index</em> can be used whenever persistence pairs are cancelled in order of increasing ratio (option <em><a href="http://localhost/dotclear/index.php?post/mse#nsig">-nsig</a></em> in <em>mse</em>), and <em>parent_index</em> whenever persistence pairs are cancelled in order of increasing difference (option <em><a href="http://localhost/dotclear/index.php?post/mse#cut">-cut</a></em> in <em>mse</em>).</li> </ul> <p><br /></p> <ul> <li>-<strong>field_value</strong> / <strong>log_field_value</strong> (Nodes and Segments) :<br /> The value of the field and it logarithm. For segments, the suffix <em>_p1</em> and <em>_p2</em> is added to indicate which extremity of the segment it corresponds to.</li> </ul> <p><br /></p> <ul> <li>-<strong>cell</strong> (Nodes and Segments) :<br /> The <em>type</em> and <em>index</em> of the cell corresponding to the node / segment in the initial network (i.e. from which the skeleton was computed). For segments, the suffix <em>_p1</em> and <em>_p2</em> 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.</li> </ul> <p><br /></p> <ul> <li>-<strong>robustness</strong> / <strong>robustness_ratio</strong> (Nodes and Segments) :<br /> The robustness of the node / segment. Robustness is a local measure of how contrasted the critical point / filament is with respect to its <em>local</em> background, and it is therefore a good way to select only highly contrasted subsets of the filamentary structures (see option <em><a href="http://localhost/dotclear/index.php?post/skelconv#trim">-trimBelow</a></em> in <em><a href="http://localhost/dotclear/index.php?post/skelconv">skelconv</a></em>). 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.</li> </ul> <p><br /></p> <ul> <li>- <strong>type</strong> (Segments only) :<br /> 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 <em>D-1</em>.</li> </ul>