Patterns in static

Apophenia

Macros | Typedefs | Functions | Variables
apop_data.c File Reference

Macros

#define Set_gsl_handler   gsl_error_handler_t *prior_handler = gsl_set_error_handler(apop_gsl_error);
 
#define Unset_gsl_handler   gsl_set_error_handler(prior_handler);
 

Typedefs

typedef int(* apop_fn_ir )(apop_data *, void *)
 

Functions

apop_dataapop_data_alloc (const size_t size1, const size_t size2, const int size3)
 
apop_dataapop_data_calloc (const size_t size1, const size_t size2, const int size3)
 
apop_dataapop_matrix_to_data (gsl_matrix *m)
 
apop_dataapop_vector_to_data (gsl_vector *v)
 
void apop_text_free (char ***freeme, int rows, int cols)
 
char apop_data_free_base (apop_data *freeme)
 
void apop_data_memcpy (apop_data *out, const apop_data *in)
 
apop_dataapop_data_copy (const apop_data *in)
 
apop_dataapop_data_stack (apop_data *m1, apop_data *m2, char posn, char inplace)
 
apop_data ** apop_data_split (apop_data *in, int splitpoint, char r_or_c)
 
void apop_data_rm_columns (apop_data *d, int *drop)
 
apop_dataapop_data_prune_columns_base (apop_data *d, char **colnames)
 
double * apop_data_ptr (apop_data *data, int row, int col, const char *rowname, const char *colname, const char *page)
 
double apop_data_get (const apop_data *data, size_t row, int col, const char *rowname, const char *colname, const char *page)
 
void apop_gsl_error_for_set (const char *reason, const char *file, int line, int gsl_errno)
 
int apop_data_set (apop_data *data, size_t row, int col, const double val, const char *colname, const char *rowname, const char *page)
 
int apop_data_set_row (apop_data *d, apop_data *row, int row_number)
 
void apop_data_add_named_elmt (apop_data *d, char *name, double val)
 
void apop_data_add_names_base (apop_data *d, const char type, char const **names)
 
int apop_text_add (apop_data *in, const size_t row, const size_t col, const char *fmt,...)
 
apop_dataapop_text_alloc (apop_data *in, const size_t row, const size_t col)
 
apop_dataapop_data_transpose (apop_data *in, char transpose_text, char inplace)
 
gsl_matrix * apop_matrix_realloc (gsl_matrix *m, size_t newheight, size_t newwidth)
 
gsl_vector * apop_vector_realloc (gsl_vector *v, size_t newheight)
 
apop_dataapop_data_get_page (const apop_data *data, const char *title, const char match)
 
apop_dataapop_data_add_page (apop_data *dataset, apop_data *newpage, const char *title)
 
apop_dataapop_data_rm_page (apop_data *data, const char *title, const char free_p)
 
apop_dataapop_data_rm_rows (apop_data *in, int *drop, apop_fn_ir do_drop, void *drop_parameter)
 

Variables

char * apop_nul_string = ""
 

Detailed Description

The apop_data structure joins together a gsl_matrix, apop_name, and a table of strings.

Function Documentation

void apop_data_add_named_elmt ( apop_data d,
char *  name,
double  val 
)

A convenience function to add a named element to a data set. Many of Apophenia's testing procedures use this to easily produce a column of named parameters. It is public as a convenience.

Parameters
dThe apop_data structure. Must not be NULL, but may be blank (as per allocation via apop_data_alloc ( ) ).
nameThe name to add
valthe value to add to the set.
  • I use the position of the last non-empty row name to know where to put the value. If there are two names in the data set, then I will put the new name in the third name slot and the data in the third slot in the vector. If you use this function from start to finish in building your list, then you'll be fine.
  • If the vector is too short (or NULL), I will call apop_vector_realloc internally to make space.
  • This fits well with the defaults for apop_data_get. An example:
1 apop_data *list = apop_data_alloc();
2 apop_data_add_named_elmt(list, "height", 165);
3 apop_data_add_named_elmt(list, "weight", 60);
4 
5 double height = apop_data_get(list, .rowname="height");
apop_data* apop_data_add_page ( apop_data dataset,
apop_data newpage,
const char *  title 
)

Add a page to a apop_data set. It gets a name so you can find it later.

Parameters
datasetThe input data set, to which a page will be added.
newpageThe page to append
titleThe name of the new page. Remember, this is truncated at 100 characters.
Returns
The new page. I post a warning if I am appending or appending to a NULL data set and apop_opts.verbose >=1 .
  • Some data is fundamentally multi-page; an optimization search over multi-page parameters would search the space given by all pages, for example. Also, pages may be appended as output or auxiliary information, such as covariances—an MLE would not search over these elements. Generally, any page with a name in XML-ish brackets, such as <Covariance>, will be considered informational and ignored by search routines, missing data routines, et cetera. This is achieved by a rule in apop_data_pack and apop_data_unpack.

Here is a toy example that establishes a baseline data set, adds a page, modifies it, and then later retrieves it.

1 apop_data *d = apop_data_alloc(10, 10, 10); //the base data set.
2 apop_data *a_new_page = apop_data_add_page(d, apop_data_alloc(2,2), "new 2 x 2 page");
3 gsl_vector_set_all(a_new_page->matrix, 3);
4 
5 //later:
6 apop_data *retrieved = apop_data_get_page(d, "new", 'r'); //use regexes, not literal match.
7 apop_data_show(retrieved); //print a 2x2 grid of 3s.
apop_data* apop_data_alloc ( const size_t  size1,
const size_t  size2,
const int  size3 
)

Allocate a apop_data structure, to be filled with data.

  • The typical case is three arguments, like apop_data_alloc(2,3,4): vector size, matrix rows, matrix cols. If the first argument is zero, you get a NULL vector.
  • Two arguments, apop_data_alloc(2,3), would allocate just a matrix, leaving the vector NULL.
  • One argument, apop_data_alloc(2), would allocate just a vector, leaving the matrix NULL.
  • Zero arguments, apop_data_alloc(), will produce a basically blank set, with out->matrix==out->vector==NULL.

For allocating the text part, see apop_text_alloc.

The weights vector is set to NULL. If you need it, allocate it via

1 d->weights = gsl_vector_alloc(row_ct);
See also
apop_data_calloc
Returns
The apop_data structure, allocated and ready.
Exceptions
out->error=='a'Allocation error. The matrix, vector, or names couldn't be malloced, which probably means that you requested a very large data set.
  • An apop_data struct, by itself, is about 72 bytes. If I can't allocate that much memory, I return NULL. But if even this much fails, your computer may be on fire and you should go put it out.
apop_data* apop_data_calloc ( const size_t  size1,
const size_t  size2,
const int  size3 
)

Allocate a apop_data structure, to be filled with data; set everything in the allocated portion to zero. See apop_data_alloc for details.

Returns
The apop_data structure, allocated and zeroed out.
Exceptions
out->error=='m'malloc error; probably out of memory.
See also
apop_data_alloc
apop_data* apop_data_copy ( const apop_data in)

Copy one apop_data structure to another. That is, all data is duplicated.

Basically a front-end for apop_data_memcpy for those who prefer this sort of syntax.

Unlike apop_data_memcpy, I do follow the more pointer.

Parameters
inthe input data
Returns
a structure that this function will allocate and fill. If input is NULL, then this will be NULL.
Exceptions
out.error='a'Allocation error.
out.error='c'Cyclic link: D->more == D (may be later in the chain, e.g., D->more->more = D->more) You'll have only a partial copy.
out.error='d'Dimension error; should never happen.
out.error='p'Missing part error; should never happen.
  • If the input data set has an error, then I will copy it anyway, including the error flag (which might be overwritten). I print a warning if the verbosity level is >=1.
char apop_data_free_base ( apop_data freeme)

Free the elements of the given apop_data set and then the apop_data set itself. Intended to be used by apop_data_free, a macro that calls this to free elements, then sets the value to NULL.

  • apop_data_free is a macro that calls this function and, on success, sets the input pointer to NULL. For typical cases, that's slightly more useful than this function.
Exceptions
freeme.error='c'Circular linking is against the rules. If freeme->more == freeme, then I set freeme.error='c' and return. If you send in a structure like A -> B -> B, then both data sets A and B will be marked.
Returns
0 on OK, 'c' on error.
apop_data* apop_data_get_page ( const apop_data data,
const char *  title,
const char  match 
)

It's good form to get a page from your data set by name, because you may not know the order for the pages, and the stepping through makes for dull code anyway (apop_data *page = dataset; while (page->more) page= page->more;).

Parameters
dataThe apop_data set to use. No default; if NULL, gives a warning if apop_opts.verbose >=1 and returns NULL.
titleThe name of the page to retrieve. Default="Info", which is the name of the page of additional estimation information returned by estimation routines (log likelihood, status, AIC, BIC, confidence intervals, ...).
matchIf 'c', case-insensitive match (via strcasecmp); if 'e', exact match, if 'r' regular expression substring search (via apop_regex). Default='c'.
Returns
The page whose title matches what you gave me. If I don't find a match, return NULL.
void apop_data_memcpy ( apop_data out,
const apop_data in 
)

Copy one apop_data structure to another.

This function does not allocate the output structure or the vector, matrix, text, or weights elements—I assume you have already done this and got the dimensions right. I will assert that there is at least enough room in the destination for your data, and fail if the copy would write more elements than there are bins.

  • If you want space allocated or are unsure about dimensions, use apop_data_copy.
  • If both in and out have a more pointer, also copy subsequent page(s).
  • You can use the subsetting macros, Apop_r or Apop_rs, to copy within a data set:
1 //Copy the contents of row i of mydata to row j.
2 apop_data *fromrow = Apop_r(mydata, i);
3 apop_data *torow = Apop_r(mydata, j);
4 apop_data_memcpy(torow, fromrow);
5 
6 // or just
7 apop_data_memcpy(Apop_r(mydata, i), Apop_r(mydata, j));
Parameters
outA structure that this function will fill. Must be preallocated with the appropriate sizes.
inThe input data.
Exceptions
out.error='d'Dimension error; couldn't copy.
out.error='p'Part missing; e.g., in->matrix exists but out->matrix doesn't; couldn't copy.
apop_data* apop_data_prune_columns_base ( apop_data d,
char **  colnames 
)

Keep only the columns of a data set that you name. This is the function called internally by the apop_data_prune_columns macro. In most cases, you'll want to use that macro. An example of the two uses demonstrating the difference:

1 apop_data_prune_columns(d, "mean", "median");
2 
3 char *list[] = {"mean", "median", NULL};
4 apop_data_prune_columns_base(d, list);
Parameters
dThe data set to prune.
colnamesA null-terminated list of names to retain (i.e. the columns that shouldn't be pruned out).
Returns
A pointer to the input data set, now pruned.
void apop_data_rm_columns ( apop_data d,
int *  drop 
)

Remove the columns set to one in the drop vector. The returned data structure looks like it was modified in place, but the data matrix and the names are duplicated before being pared down, so if your data is taking up more than half of your memory, this may not work.

Parameters
dthe apop_data structure to be pared down.
dropan array of ints. If use[7]==1, then column seven will be cut from the output. A reminder: calloc(in->size2 , sizeof(int)) will fill your array with zeros on allocation, and memset(use, 1, in->size2 * sizeof(int)) will quickly fill an array of ints with nonzero values.
apop_data* apop_data_rm_page ( apop_data data,
const char *  title,
const char  free_p 
)

Remove the first page from an apop_data set that matches a given name.

Parameters
dataThe input data set, to which a page will be added. No default. If NULL, I return silently if apop_opts.verbose < 1 ; print an error otherwise.
titleThe case-insensitive name of the page to remove. Default: "Info"
free_pIf 'y', then apop_data_free the page. Default: 'y'.
Returns
If not freed, a pointer to the apop_data page that I just pulled out. Thus, you can use this to pull a single page from a data set. I set that page's more pointer to NULL, to minimize any confusion about more-than-linear linked list topologies. If free_p=='y' (the default) or the page is not found, return NULL.
  • I don't check the first page, so there's no concern that the head of your list of pages will move. Again, the intent of the ->more pointer in the apop_data set is not to fully implement a linked list, but primarily to allow you to staple auxiliary information to a main data set.
  • If I don't find the page you want, I return NULL, and print a message if apop_opts.verbose >= 1.
apop_data* apop_data_rm_rows ( apop_data in,
int *  drop,
apop_fn_ir  do_drop,
void *  drop_parameter 
)

Remove the rows set to one in the drop vector or for which the do_drop function returns one.

Parameters
inthe apop_data structure to be pared down
dropa vector with as many elements as the max of the vector, matrix, or text parts of in, with a one marking those columns to be removed.
do_dropA function that returns one for rows to drop and zero for rows to not drop. A sample function:
1 int your_drop_function(apop_data *onerow, void *extra_param){
2  return gsl_isnan(apop_data_get(onerow)) || !strcmp(onerow->text[0][0], "Uninteresting data point");
3 }
apop_data_rm_rows uses Apop_r to get a subview of the input data set of height one (and since all the default arguments default to zero, you don't have to write out things like apop_data_get (onerow, .row=0, .col=0), which can help to keep things readable).
drop_parameterIf your do_drop function requires additional input, put it here and it will be passed through.
Returns
Returns a pointer to the input data set, now pruned.
  • If all the rows are to be removed, then you will wind up with the same apop_data set, with NULL vector, matrix, weight, and text. Therefore, you may wish to check for NULL elements after use. I remove rownames, but leave the other names, in case you want to add new data rows.
  • The typical use is to provide only a list or only a function. If both are NULL, I return without doing anything, and print a warning if apop_opts.verbose >=1. If you provide both, I will drop the row if either the vector has a one in that row's position, or if the function returns a nonzero value.
  • This function uses the Designated initializers syntax for inputs.
int apop_data_set_row ( apop_data d,
apop_data row,
int  row_number 
)

Now that you've used Apop_r to pull a row from an apop_data set, this function lets you write that row to another position in the same data set or a different data set entirely.

The set written to must have the same form as the original:

  • a vector element has to be present if one existed in the original,
  • same for the weights vector,
  • the matrix in the destination has to have as many columns as in the original, and
  • the text has to have a row long enough to hold the original
  • If the row to be written to already has a rowname, it is overwritten. If d->names->rowct == row_number (all rows up to row_number have row names), then extend the list of row names by one to add the new name. Else, don't add the row name.
  • Column names (of all types) aren't touched. Maybe use apop_data_copy or apop_name_copy if you need to copy these names.

If any of the source elements are NULL, I won't bother to check that element in the destination.

Returns
0=OK, -1=error (probably a source/destination size mismatch).
  • The error codes for out-of-bounds errors are thread-safe iff you are have a C11-compliant compiler (thanks to the _Thread_local keyword) or a version of GCC with the __thread extension enabled.
apop_data** apop_data_split ( apop_data in,
int  splitpoint,
char  r_or_c 
)

Split one input apop_data structure into two.

For the opposite operation, see apop_data_stack.

Parameters
inThe apop_data structure to split
splitpointThe index of what will be the first row/column of the second data set. E.g., if this is -1 and r_or_c=='c', then the whole data set will be in the second data set; if this is the length of the matrix then the whole data set will be in the first data set. Another way to put it is that splitpoint will equal the number of rows/columns in the first matrix (unless it is -1, in which case the first matrix will have zero rows, or it is greater than the matrix's size, in which case it will have as many rows as the original).
r_or_cIf this is 'r' or 'R', then put some rows in the first data set and some in the second; of 'c' or 'C', split columns into first and second data sets.
Returns
An array of two apop_data sets. If one is empty then a NULL pointer will be returned in that position. For example, for a data set of 50 rows, apop_data **out = apop_data_split(data, 100, 'r') sets out[0] = apop_data_copy(data) and out[1] = NULL.
  • When splitting at a row, the text is also split.
  • more pointer is ignored.
  • The apop_data->vector is taken to be the -1st element of the matrix.
  • Weights will be preserved. If splitting by rows, then the top and bottom parts of the weights vector will be assigned to the top and bottom parts of the main data set. If splitting by columns, identical copies of the weights vector will be assigned to both parts.
  • Data is copied, so you may want to call apop_data_free(in) after this.
apop_data* apop_data_stack ( apop_data m1,
apop_data m2,
char  posn,
char  inplace 
)

Put the first data set either on top of or to the left of the second data set.

The fn returns a new data set, meaning that at the end of this function, until you apop_data_free() the original data sets, you will be taking up twice as much memory. Plan accordingly.

For the opposite operation, see apop_data_split.

Parameters
m1the upper/rightmost data set (default = NULL)
m2the second data set (default = NULL)
posnIf 'r', stack rows of m1's matrix above rows of m2's
if 'c', stack columns of m1's matrix to left of m2's
(default = 'r')
inplaceIf 'y', use apop_matrix_realloc and apop_vector_realloc to modify m1 in place; see the caveats on those function. Otherwise, allocate a new vector, leaving m1 unmolested. (default='n')
Returns
The stacked data, either in a new apop_data set or m1
Exceptions
out->error=='a'Allocation error.
out->error=='d'Dimension error; couldn't make a complete copy.
  • If m1 or m2 are NULL, this returns a copy of the other element, and if both are NULL, you get NULL back (except if m2 is NULL and inplace is 'y', where you'll get the original m1 pointer back)
  • Text is handled as you'd expect: If 'r', one set of text is stacked on top of the other [number of columns must match]; if 'c', one set of text is set next to the other [number of rows must match].
  • more is ignored.
  • If stacking rows on rows, the output vector is the input vectors stacked accordingly. If stacking columns by columns, the output vector is just a copy of the vector of m1 and m2->vector doesn't appear in the output at all.
  • The same rules for dealing with the vector(s) hold for the vector(s) of weights.
  • Names are a copy of the names for m1, with the names for m2 appended to the row or column list, as appropriate.
  • This function uses the Designated initializers syntax for inputs.
apop_data* apop_data_transpose ( apop_data in,
char  transpose_text,
char  inplace 
)

Transpose the matrix and text elements of the input data set, including the row/column names.

The vector and weights elements of the input data set are completely ignored (but see also apop_vector_to_matrix, which can convert a vector to a 1 X N matrix.) If copying, these other elements won't be present; if .inplace='y', it is up to you to handle these not-transposed elements correctly.

Parameters
inThe input apop_data set. If NULL, I return NULL. Default is NULL.
transpose_textIf 'y', then also transpose the text element. Default is 'y'.
inplaceIf 'y', transpose the input in place; if 'n', produce a transposed copy, leaving the original untouched. Due to how gsl_matrix_transpose_memcpy works, a copy will still be made, then copied to the original location. Default is 'y'.
Returns
If inplace=='n', a newly alloced apop_data set, with the appropriately transposed matrix and/or text. The vector and weights elements will be NULL. If transpose_text='n', then the text element of the output set will also be NULL.
if inplace=='y', a pointer to the original data set, with matrix and (if transpose_text='y') text transposed and vector and weights left in place untouched.
  • Row names are written to column names of the output matrix, text, or both (whichever is not empty in the input).
  • If only the matrix or only the text have names, then the one set of names is written to the row names of the output.
  • If both matrix column names and text column names are present, text column names are lost.
  • if you have a gsl_matrix with no names or text, you may prefer to use gsl_matrix_transpose_memcpy.
gsl_matrix* apop_matrix_realloc ( gsl_matrix *  m,
size_t  newheight,
size_t  newwidth 
)

This function will resize a gsl_matrix to a new height or width.

Data in the matrix will be retained. If the new height or width is smaller than the old, then data in the later rows/columns will be cropped away (in a non–memory-leaking manner). If the new height or width is larger than the old, then new cells will be filled with garbage; it is your responsibility to zero out or otherwise fill new rows/columns before use.

Warning I: Using this function is basically bad form—especially when used in a for loop that adds a column each time. A large number of reallocs can take a noticeable amount of time. You are thus encouraged to make an effort to determine the size of your data beforehand.

Warning II: The gsl_matrix is a versatile struct that can represent submatrices and other cuts from parent data. I can't deal with those, and check for such situations beforehand. [Besides, resizing a portion of a parent matrix makes no sense.]

Parameters
mThe already-allocated matrix to resize. If you give me NULL, this becomes equivalent to gsl_matrix_alloc
newheight,newwidthThe height and width you'd like the matrix to be.
Returns
m, now resized
apop_data* apop_matrix_to_data ( gsl_matrix *  m)

Deprecated; please do not use. Just use a compound literal:

1 //Given:
2 gsl_vector *v;
3 gsl_matrix *m;
4 
5 // Then this form wraps the elements into \ref apop_data structs. Note that
6 // these are not pointers: they're automatically allocated and therefore
7 // the extra memory use for the wrapper is cleaned up on exit from scope.
8 
9 apop_data *dv = &(apop_data){.vector=v};
10 apop_data *dm = &(apop_data){.matrix=m};
11 
12 apop_data *v_dot_m = apop_dot(dv, dm);
13 
14 //Here is a macro to hide C's ugliness:
15 #define As_data(...) (&(apop_data){__VA_ARGS__})
16 
17 apop_data *v_dot_m2 = apop_dot(As_data(.vector=v), As_data(.matrix=m));
18 
19 //The wrapped object is an automatically-allocated structure pointing to the
20 //original data. If it needs to persist or be separate from the original,
21 //make a copy:
22 apop_data *dm_copy = apop_data_copy(As_data(.vector=v, .matrix=m));
int apop_text_add ( apop_data in,
const size_t  row,
const size_t  col,
const char *  fmt,
  ... 
)

Add a string to the text element of an apop_data set. If you send me a NULL string, I will write the value of apop_opts.nan_string in the given slot. If there is already something in that slot, that string is freed, preventing memory leaks.

Parameters
inThe apop_data set, that already has an allocated text element.
rowThe row
colThe col
fmtThe text to write.
...You can use a printf-style fmt and follow it with the usual variables to fill in.
Returns
0=OK, -1=error (probably out-of-bounds)
  • UTF-8 or ASCII text is correctly handled.
  • Apophenia follows a general rule of not reallocating behind your back: if your text matrix is currently of size (3,3) and you try to put an item in slot (4,4), then I display an error rather than reallocating the text matrix.
  • Resizing a text matrix is annoying in C, so note that apop_text_alloc will reallocate to a new size if you need. For example, this code will fill the diagonals of the text array with a message, resizing as it goes:
  • The string added is a copy (via asprintf), not a pointer to the input(s).
  • If there had been a string at the grid point you are writing to, the old one is effectively lost when the new one is placed. So, I free the old string to prevent leaks if necessary. Remember this if you had other pointers aliasing that string, in which case you may as well avoid this function and just use asprintf(&(your_dataset->text[row][col]), "your string").
1 apop_data *list = (something already allocated.);
2 for (int n=0; n < 10; n++){
3  apop_text_alloc(list, n+1, n+1);
4  apop_text_add(list, n, n, "This is cell (%i, %i)", n, n);
5 }
apop_data* apop_text_alloc ( apop_data in,
const size_t  row,
const size_t  col 
)

This allocates an array of strings and puts it in the text element of an apop_data set.

If the text element already exists, then this is effectively a realloc function, reshaping to the size you specify.

Parameters
inAn apop_data set. It's OK to send in NULL, in which case an apop_data set with NULL matrix and vector elements is returned.
rowthe number of rows of text.
colthe number of columns of text.
Returns
A pointer to the relevant apop_data set. If the input was not NULL, then this is a repeat of the input pointer.
Exceptions
out->error=='a'Allocation error.
void apop_text_free ( char ***  freeme,
int  rows,
int  cols 
)

Free a matrix of chars* (i.e., a char***). This is the form of the text element of the apop_data set, so you can use this for:

1 apop_text_free(yourdata->text, yourdata->textsize[0], yourdata->textsize[1]);

This is what apop_data_free uses internally.

gsl_vector* apop_vector_realloc ( gsl_vector *  v,
size_t  newheight 
)

This function will resize a gsl_vector to a new length.

Data in the vector will be retained. If the new height is smaller than the old, then data at the end of the vector will be cropped away (in a non–memory-leaking manner). If the new height is larger than the old, then new cells will be filled with garbage; it is your responsibility to zero out or otherwise fill them before use.

Warning I: Using this function is basically bad form—especially when used in a for loop that adds an element each time. A large number of reallocs can take a noticeable amount of time. You are thus encouraged to make an effort to determine the size of your data beforehand.

Warning II: The gsl_vector is a versatile struct that can represent subvectors, matrix columns and other cuts from parent data. I can't deal with those, and check for such situations beforehand. [Besides, resizing a portion of a parent matrix makes no sense.]

Parameters
vThe already-allocated vector to resize. If you give me NULL, this is equivalent to gsl_vector_alloc
newheightThe height you'd like the vector to be.
Returns
v, now resized
apop_data* apop_vector_to_data ( gsl_vector *  v)

Deprecated; please do not use. Just use a compound literal, as in the code sample in the documentation for apop_matrix_to_data.

Autogenerated by doxygen on Wed Oct 15 2014 (Debian 0.999b+ds3-2).