![]() |
|
Go to the source code of this file.
| typedef struct _apop_model apop_model |
A description of a parametrized statistical model, including the input settings and the output parameters, expected values, et cetera. The full declaration is given in the _apop_model page, see the longer discussion on the models page, or see the apop_ols page for a sample program that uses an apop_model.
| void apop_data_add_named_elmt | ( | apop_data * | d, | |
| char * | name, | |||
| double | val | |||
| ) |
A convenience function to add a named element to a data vector. For example, many of the testing procedures use this to easily produce a column of named parameters.
| d | The apop_data structure. Must not be NULL. | |
| name | The name to add | |
| val | the value to add to the vector. |
NULL), I will call apop_vector_realloc internally to make space in the vector. Split one input apop_data structure into two.
For the opposite operation, see apop_data_stack.
The apop_data->vector is taken to be the -1st element of the matrix.
| in | The apop_data structure to split | |
| splitpoint | The 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_c | If this is 'r' or 'R', then cleave the rows; of 'c' or 'C' cleave the columns. |
apop_data **out = apop_data_split(data, 100, 'r') sets out[0] = apop_data_copy(data) and out[1] = NULL.Transpose the matrix element of the input apop_data set, including the row/column names. The vector and text elements of the input data set are completely ignored.
This is really just a friendly wrapper for gsl_matrix_transpose_memcpy; if you have a gsl_matrix with no names, you may prefer to just use that function.
| in | The input apop_data set. |
NULL. | 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 repsonsibility 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.]
| m | The already-allocated matrix to resize. If you give me NULL, this becomes equivalent to gsl_matrix_alloc | |
| newheight,newwidth | The height and width you'd like the matrix to be. |
| apop_data* apop_matrix_to_data | ( | gsl_matrix * | m | ) |
Wrap an apop_data structure around an existing gsl_matrix. The matrix is not copied, but is pointed to by the new apop_data struct.
| m | The existing matrix you'd like to turn into an apop_data structure. |
matrix pointer points to the input matrix. The rest of the struct is basically blank. | int apop_name_add | ( | apop_name * | n, | |
| char * | add_me, | |||
| char | type | |||
| ) |
Adds a name to the apop_name structure. Puts it at the end of the given list.
| n | An existing, allocated apop_name structure. | |
| add_me | A string. If NULL, do nothing; return -1. | |
| type | 'r': add a row name 'c': add a column name 't': add a text category name 'h': add a title (or a header. 't' is taken). 'v': add (or overwrite) the vector name |
| apop_name* apop_name_alloc | ( | void | ) |
Allocates a name structure
Copy one apop_name structure to another. That is, all data is duplicated. Usage:
apop_name *out = apop_name_copy(in);
| in | the input names |
Deprecated. Use apop_name_stack.
| int apop_name_find | ( | apop_name * | n, | |
| char * | in, | |||
| char | type | |||
| ) |
Finds the position of an element in a list of names.
The function uses case-insensitive regular expressions to search.
For example, "p.val.*" will match "P value", "p.value", and "p values".
| n | the apop_name object to search. | |
| in | the name you seek; see above. | |
| type | 'c', 'r', or 't'. Default is 'c'. |
findme. If 'c', then this may be -1, meaning the vector name. If not found, returns -2. | void apop_name_print | ( | apop_name * | n | ) |
Prints the given list of names to STDOUT
| n | the apop_name structure |
Append one list of names to another.
Notice that if the first list is NULL, then this is a copy function. If the second is NULL, it is a no-op.
| n1 | The first set of names (no default, must not be NULL) | |
| nadd | The second set of names, which will be appended after the first. (no default, if NULL, a no-op) | |
| type1 | Either 'c', 'r', 't', or 'v' stating whether you are merging the columns, rows, or text. If 'v', then ignore typeadd and just overwrite the target vector name with the source name. (default = 'r') | |
| typeadd | Either 'c', 'r', 't', or 'v' stating whether you are merging the columns, rows, or text. If 'v', then overwrite the target with the source vector name. (default = type1) |
| void 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 string "NaN" in the given slot.
| in | The apop_data set, that already has an allocated text element. | |
| row | The row | |
| col | The col | |
| fmt | The text to write. | |
| ... | You can use a printf-style fmt and follow it with the usual variables to fill in. |
This allocates an array of strings and puts it in the text element of an apop_data set.
| in | An 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. | |
| row | the number of rows of text. | |
| col | the number of columns of text. |
NULL, then this is a repeat of the input pointer. | 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:
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 in the bottom 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 repsonsibility 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.]
| v | The already-allocated vector to resize. If you give me NULL, this is equivalent to gsl_vector_alloc | |
| newheight | The height you'd like the vector to be. |
| apop_data* apop_vector_to_data | ( | gsl_vector * | v | ) |
Here are where the options are initially set.