Patterns in static

Apophenia

Printing to the screen or a text file

Modules

Functions


Detailed Description

Most functions print only to the screen, but the matrix and vector printing functions will let you print to a text file as well. The presumption is that statistic estimates are for your own consumption, while you are printing a matrix for import into another program.


Function Documentation

void apop_data_show ( const apop_data in  ) 

This function prettyprints the apop_data set to a screen.

This takes a lot of machinery. I write every last element to a text array, then measure column widths, then print to screen with padding to guarantee that everything lines up. There's no way to have the first element of a column line up with the last unless you interrogate the width of every element in the column, so printing columns really can't be a one-pass process.

So, I produce an apop_data set with no numeric elements and a text element to be filled with the input data set, and then print that. That means that I'll be using (more than) twice the memory to print this. If this is a problem, you can use Assorted printing functions to dump your data to a text file, and view the text file, or print subsets.

For more machine-readable printing, see Assorted printing functions.

apop_data* apop_data_summarize ( apop_data indata  ) 

Put summary information about the columns of a table (mean, std dev, variance) in a table.

Parameters:
indata The table to be summarized. An apop_data structure.
Returns:
An apop_data structure with one row for each column in the original table, and a column for each summary statistic. May have a weights element.
Todo:

At the moment, only gives the mean, standard deviation, and variance of the data in each column; should give more in the near future.

We should probably let this summarize rows as well.

void apop_histogram_plot ( apop_model hist,
Output_declares   
)

This function will plot a histogram model. Compare with apop_plot_histogram, which is a convenience function to plot a gsl_vector as a histogram.

The function respects the output_type option, so code like:

apop_opts.output_type = 'p';
apop_opts.output_pipe = popen("/usr/bin/gnuplot", "w");
apop_model *m = apop_estimate(data, apop_histogram);
apop_plot_histogram(m);

will print directly to Gnuplot.

Parameters:
hist A parametrized apop_model holding the histogram. (No default. Must not be NULL.)

This function uses the Designated initializers syntax for inputs.

void apop_model_show ( apop_model print_me  ) 

Print the results of an estimation. If your model has a print method, then I'll use that, else I'll use a default setup.

Your print method can use both by masking itself for a second:

void print_method(apop_model *in){
  void *temp = in->print;
  in->print = NULL;
  apop_model_show(in);
  in->print = temp;

  printf("Additional info:\n");
  ...
}
void apop_plot_histogram ( gsl_vector *  data,
size_t  bins,
Output_declares   
)

This function will take in a gsl_vector of data and put out a histogram. Compare with apop_histogram_plot, which plots an estimated apop_histogram model.

The function respects the output_type option, so code like:

f   = popen("/usr/bin/gnuplot", "w");
apop_opts.output_type = 'p';
apop_opts.output_pipe = f;
apop_plot_histogram(data, 100, NULL);

will print directly to Gnuplot.

Parameters:
data A gsl_vector holding the data. Do not pre-sort or bin; this function does that for you. (no default, must not be NULL)
bins The number of bins in the output histogram (default = MAX(10, data->size/20); denominator subject to future adjustment)

This function uses the Designated initializers syntax for inputs.

void apop_plot_lattice ( const apop_data d,
Output_declares   
)

This produces a Gnuplot file that will produce an array of 2-D plots, one for each pair of columns in the data set. Along the diagonal is a plot of the variable against itself---a density plot of the variable.

Parameters:
d The data set whose (matrix) columns will be compared. (No default, must not be NULL.)
lattice.png

A lattice showing three variables graphed against each other.

  • This function uses the Designated initializers syntax for inputs.
  • This function respects the global apop_opts.output_type and apop_opts.output_pipe variables; see the legible output section of the outline for details.
void apop_plot_line_and_scatter ( apop_data data,
apop_model est,
Output_declares   
)

Prep for gnuplot one of those cute scatterplots with a regression line through it.

Currently, you only get two dimensions.

Set the global apop_opts.output_name to the filename you want before running this. It appends instead of overwriting, so you can prep the file if you want; see sample code. [to overwrite a file, just remove it first with the standard C function remove("filename");]

Parameters:
data This is a copy of what you'd sent to the regression fn. That is, the first column is the dependent variable and the second is the independent. That is, what will be on Y axis is the first column, and what is on the X axis is the second. Custom for regressions and custom for graphs just clash on this one.
est The apop_model structure your regression function gave you. (if NULL, I'll estimate an OLS model for you).

The sample program below will pull data from a database (ridership at the Silver Spring, MD Metro station; get the database in the {Modeling with Data} sample code, at http://modelingwithdata.org/appendices.html), then runs OLS, and produce a gnuplot file to write to a file named "scatter.eps". You can run the result through gnuplot via gnuplot scatter.gplot, and if you don't like the results, you have the gnuplot file ("scatter.gplot") on hand for modifications.

#include <apop.h>

int main(){
apop_data   *data, *data_copy;
apop_model  *est;
FILE        *f;
char        outfile[]   = "scatter.gplot";

    apop_db_open("data-metro.db");
    data = apop_query_to_data("select riders, year from riders where station like 'Silver%%'");
    apop_db_close();

    //The regression destroys your data, so copy it first.
    data_copy   = apop_data_copy(data);

    //Run OLS, display results on terminal
    est  = apop_estimate(data, apop_OLS);
    apop_model_show(est);

    //Prep the file with a header, then call the function.
    f    = fopen(outfile, "w");
    fprintf(f,"set term postscript;\n set output \"scatter.eps\"\n set yrange [0:*]\n");
    apop_plot_line_and_scatter(data_copy, est, .output_pipe=f);
    fclose(f);

    return 0;
}
  • This function uses the Designated initializers syntax for inputs.
  • This function respects the global apop_opts.output_type and apop_opts.output_pipe variables; see the legible output section of the outline for details.

SourceForge.net Logo

Autogenerated by doxygen on 23 Nov 2009.