-
Notifications
You must be signed in to change notification settings - Fork 31
/
apop.m4.h
1704 lines (1388 loc) · 85.7 KB
/
apop.m4.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/** \file */
/* Copyright (c) 2005--2014 by Ben Klemens. Licensed under the GPLv2; see COPYING. */
/* Here are the headers for all of apophenia's functions, typedefs, static variables and
macros. All of these begin with the apop_ (or Apop_ or APOP_) prefix.
There used to be a series of sub-headers, but they never provided any serious
benefit. Please use your text editor's word-search feature to find any elements you
may be looking for. About a third of the file is comments and doxygen documentation,
so syntax highlighting that distinguishes code from comments will also help to make
this more navigable.*/
/** \defgroup all_public Public functions, structs, and types
\addtogroup all_public
@{
*/
#pragma once
#ifdef __cplusplus
extern "C" {
#endif
/** \cond doxy_ignore */
#ifndef _GNU_SOURCE
#define _GNU_SOURCE //for asprintf
#endif
#include <assert.h>
#include <signal.h> //raise(SIGTRAP)
#include <string.h>
#include <gsl/gsl_rng.h>
#include <gsl/gsl_matrix.h>
//////Optional arguments
/* A means of providing more script-like means of sending arguments to a function.
These macros are intended as internal. If you are interested in using this mechanism
in out-of-Apophenia work, grep docs/documentation.h for optionaldetails to find notes
on how these are used (Doxygen doesn't use that page),
*/
#define apop_varad_head(type, name) type variadic_##name(variadic_type_##name varad_in)
#define apop_varad_declare(type, name, ...) \
typedef struct { \
__VA_ARGS__ ; \
} variadic_type_##name; \
apop_varad_head(type, name);
#define apop_varad_var(name, value) name = varad_in.name ? varad_in.name : (value);
#define apop_varad_link(name,...) variadic_##name((variadic_type_##name) {__VA_ARGS__})
/** \endcond */ //End of Doxygen ignore.
//////The types and functions that act on them
/** This structure holds the names of the components of the \ref apop_data set. You may never have to worry about it directly, because most operations on \ref apop_data sets will take care of the names for you.
*/
typedef struct{
char *title;
char * vector;
char ** col;
char ** row;
char ** text;
int colct, rowct, textct;
unsigned long *colhash, *rowhash, *texthash;
} apop_name;
/** The \ref apop_data structure represents a data set. See \ref dataoverview.*/
typedef struct apop_data{
gsl_vector *vector;
gsl_matrix *matrix;
apop_name *names;
char ***text;
size_t textsize[2];
gsl_vector *weights;
struct apop_data *more;
char error;
} apop_data;
/* Settings groups. For internal use only; see apop_settings.c and
settings.h for related machinery. */
typedef struct {
char name[101];
unsigned long name_hash;
void *setting_group;
void *copy;
void *free;
} apop_settings_type;
/** A statistical model. See \ref modelsec for details. */
typedef struct apop_model apop_model;
/** The elements of the \ref apop_model type, representing a statistical model. See \ref
modelsec and \ref modeldetails for use and details. */
struct apop_model{
char name[101];
int vsize, msize1, msize2, dsize;
apop_data *data;
apop_data *parameters;
apop_data *info;
void (*estimate)(apop_data * data, apop_model *params);
long double (*p)(apop_data *d, apop_model *params);
long double (*log_likelihood)(apop_data *d, apop_model *params);
long double (*cdf)(apop_data *d, apop_model *params);
long double (*constraint)(apop_data *data, apop_model *params);
int (*draw)(double *out, gsl_rng* r, apop_model *params);
void (*prep)(apop_data *data, apop_model *params);
apop_settings_type *settings;
void *more;
size_t more_size;
char error;
};
/** The global options. */
typedef struct{
int verbose; /**< Set this to zero for silent mode, one for errors and warnings. default = 0. */
char stop_on_warning; /**< See \ref debugging . */
char output_delimiter[100]; /**< The separator between elements of output tables. The default is "\t", but
for LaTeX, use "&\t", or use "|" to get pipe-delimited output. */
char input_delimiters[100]; /**< Deprecated. Please use per-function inputs to \ref apop_text_to_db and \ref apop_text_to_data. Default = "|,\t" */
char *db_name_column; /**< If not NULL or <tt>""</tt>, the name of the column in your tables that holds row names.*/
char *nan_string; /**< The string used to indicate NaN. Default: <tt>"NaN</tt>. Comparisons are case-insensitive.*/
char db_engine; /**< If this is 'm', use mySQL, else use SQLite. */
char db_user[101]; /**< Username for database login. Max 100 chars. */
char db_pass[101]; /**< Password for database login. Max 100 chars. */
FILE *log_file; /**< The file handle for the log. Defaults to \c stderr, but change it with, e.g.,
<tt>apop_opts.log_file = fopen("outlog", "w");</tt> */
#define Autoconf_no_atomics @Autoconf_no_atomics@
#if __STDC_VERSION__ > 201100L && !defined(__STDC_NO_ATOMICS__) && Autoconf_no_atomics==0
_Atomic(int) rng_seed;
#else
int rng_seed;
#endif
float version;
} apop_opts_type;
extern apop_opts_type apop_opts;
apop_name * apop_name_alloc(void);
int apop_name_add(apop_name * n, char const *add_me, char type);
void apop_name_free(apop_name * free_me);
void apop_name_print(apop_name * n);
Apop_var_declare( void apop_name_stack(apop_name * n1, apop_name *nadd, char type1, char typeadd) )
apop_name * apop_name_copy(apop_name *in);
int apop_name_find(const apop_name *n, const char *findme, const char type);
void apop_data_add_names_base(apop_data *d, const char type, char const ** names);
/** Add a list of names to a data set.
\li Use this with a list of names that you type in yourself, like
\code
apop_data_add_names(mydata, 'c', "age", "sex", "height");
\endcode
Notice the lack of curly braces around the list.
\li You may have an array of names, probably autogenerated, that you would like to
add. In this case, make certain that the last element of the array is \c NULL, and
call the base function:
\code
char **[] colnames = {"age", "sex", "height", NULL};
apop_data_add_names_base(mydata, 'c', colnames);
\endcode
But if you forget the \c NULL marker, this has good odds of segfaulting. You may prefer to use a \c for loop that inserts each name in turn using \ref apop_name_add.
\see \ref apop_name_add, although \ref apop_data_add_names will be more useful in most cases.
*/
#define apop_data_add_names(dataset, type, ...) apop_data_add_names_base((dataset), (type), (char const*[]) {__VA_ARGS__, NULL})
/** Free an \ref apop_data structure.
\li As with \c free(), it is safe to send in a \c NULL pointer (in which case the function does nothing).
\li If the \c more pointer is not \c NULL, I will free the pointed-to data set first.
If you don't want to free data sets down the chain, set <tt>more=NULL</tt> before calling this.
\li This is actually a macro (that calls \ref apop_data_free_base). It
sets \c freeme to \c NULL when it's done, because there's nothing safe you can do with the
freed location, and you can later safely test conditions like <tt>if (data) ...</tt>.
*/
#define apop_data_free(freeme) (apop_data_free_base(freeme) ? 0 : ((freeme)= NULL))
char apop_data_free_base(apop_data *freeme);
Apop_var_declare( apop_data * apop_data_alloc(const size_t size1, const size_t size2, const int size3) )
Apop_var_declare( apop_data * apop_data_calloc(const size_t size1, const size_t size2, const int size3) )
Apop_var_declare( apop_data * apop_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);
apop_data * apop_data_copy(const apop_data *in);
void apop_data_rm_columns(apop_data *d, int *drop);
void apop_data_memcpy(apop_data *out, const apop_data *in);
Apop_var_declare( double * apop_data_ptr(apop_data *data, int row, int col, const char *rowname, const char *colname, const char *page) )
Apop_var_declare( double apop_data_get(const apop_data *data, size_t row, int col, const char *rowname, const char *colname, const char *page) )
Apop_var_declare( int apop_data_set(apop_data *data, size_t row, int col, const double val, const char *rowname, const char * colname, const char *page) )
void apop_data_add_named_elmt(apop_data *d, char *name, double val);
int apop_text_set(apop_data *in, const size_t row, const size_t col, const char *fmt, ...);
apop_data * apop_text_alloc(apop_data *in, const size_t row, const size_t col);
void apop_text_free(char ***freeme, int rows, int cols);
Apop_var_declare( apop_data * apop_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);
#define apop_data_prune_columns(in, ...) apop_data_prune_columns_base((in), (char *[]) {__VA_ARGS__, NULL})
apop_data* apop_data_prune_columns_base(apop_data *d, char **colnames);
Apop_var_declare( apop_data * apop_data_get_page(const apop_data * data, const char * title, const char match) )
apop_data * apop_data_add_page(apop_data * dataset, apop_data *newpage,const char *title);
Apop_var_declare( apop_data* apop_data_rm_page(apop_data * data, const char *title, const char free_p) )
Apop_var_declare( apop_data * apop_data_rm_rows(apop_data *in, int *drop, int (*do_drop)(apop_data* ! void*), void* drop_parameter) )
//in apop_asst.c:
Apop_var_declare( apop_data * apop_model_draws(apop_model *model, int count, apop_data *draws) )
/* Convenience functions to convert among vectors (gsl_vector), matrices (gsl_matrix),
arrays (double **), and database tables */
//From vector
gsl_vector *apop_vector_copy(const gsl_vector *in);
Apop_var_declare( gsl_matrix * apop_vector_to_matrix(const gsl_vector *in, char row_col) )
//From matrix
gsl_matrix *apop_matrix_copy(const gsl_matrix *in);
Apop_var_declare( apop_data *apop_db_to_crosstab(char const*tabname, char const*row, char const*col, char const*data, char is_aggregate) )
//From array
Apop_var_declare( gsl_vector * apop_array_to_vector(double *in, int size) )
/** \cond doxy_ignore */ //Deprecated
#define apop_text_add apop_text_set
#define apop_line_to_vector apop_array_to_vector
/** \endcond */
//From text
Apop_var_declare( apop_data * apop_text_to_data(char const *text_file, int has_row_names, int has_col_names, int const *field_ends, char const *delimiters) )
Apop_var_declare( int apop_text_to_db(char const *text_file, char *tabname, int has_row_names, int has_col_names, char **field_names, int const *field_ends, apop_data *field_params, char *table_params, char const *delimiters, char if_table_exists) )
//rank data
apop_data *apop_data_rank_expand (apop_data *in);
Apop_var_declare( apop_data *apop_data_rank_compress (apop_data *in, int min_bins) )
//From crosstabs
void apop_crosstab_to_db(apop_data *in, char *tabname, char *row_col_name,
char *col_col_name, char *data_col_name);
//packing data into a vector
Apop_var_declare( gsl_vector * apop_data_pack(const apop_data *in, gsl_vector *out, char more_pages, char use_info_pages) )
Apop_var_declare( void apop_data_unpack(const gsl_vector *in, apop_data *d, char use_info_pages) )
#define apop_vector_fill(avfin, ...) apop_vector_fill_base((avfin), (double []) {__VA_ARGS__})
#define apop_data_fill(adfin, ...) apop_data_fill_base((adfin), (double []) {__VA_ARGS__})
#define apop_text_fill(dataset, ...) apop_text_fill_base((dataset), (char* []) {__VA_ARGS__, NULL})
#define apop_data_falloc(sizes, ...) apop_data_fill(apop_data_alloc sizes, __VA_ARGS__)
apop_data *apop_data_fill_base(apop_data *in, double []);
gsl_vector *apop_vector_fill_base(gsl_vector *in, double []);
apop_data *apop_text_fill_base(apop_data *data, char* text[]);
//// Models and model support functions
extern apop_model *apop_beta;
extern apop_model *apop_bernoulli;
extern apop_model *apop_binomial;
extern apop_model *apop_chi_squared;
extern apop_model *apop_dirichlet;
extern apop_model *apop_exponential;
extern apop_model *apop_f_distribution;
extern apop_model *apop_gamma;
extern apop_model *apop_improper_uniform;
extern apop_model *apop_iv;
extern apop_model *apop_kernel_density;
extern apop_model *apop_loess;
extern apop_model *apop_logit;
extern apop_model *apop_lognormal;
extern apop_model *apop_multinomial;
extern apop_model *apop_multivariate_normal;
extern apop_model *apop_normal;
extern apop_model *apop_ols;
extern apop_model *apop_pmf;
extern apop_model *apop_poisson;
extern apop_model *apop_probit;
extern apop_model *apop_t_distribution;
extern apop_model *apop_uniform;
//extern apop_model *apop_wishart;
extern apop_model *apop_wls;
extern apop_model *apop_yule;
extern apop_model *apop_zipf;
//model transformations
extern apop_model *apop_coordinate_transform;
extern apop_model *apop_composition;
extern apop_model *apop_dconstrain;
extern apop_model *apop_mixture;
extern apop_model *apop_cross;
/** Alias for the \ref apop_normal distribution, qv. */
#define apop_gaussian apop_normal
#define apop_OLS apop_ols
#define apop_PMF apop_pmf
#define apop_F_distribution apop_f_distribution
#define apop_IV apop_iv
void apop_model_free (apop_model * free_me);
Apop_var_declare( void apop_model_print (apop_model * model, FILE *output_pipe) )
void apop_model_show (apop_model * print_me); //deprecated
apop_model * apop_model_copy(apop_model *in); //in apop_model.c
apop_model * apop_model_clear(apop_data * data, apop_model *model);
apop_model * apop_estimate(apop_data *d, apop_model *m);
void apop_score(apop_data *d, gsl_vector *out, apop_model *m);
double apop_log_likelihood(apop_data *d, apop_model *m);
double apop_p(apop_data *d, apop_model *m);
double apop_cdf(apop_data *d, apop_model *m);
int apop_draw(double *out, gsl_rng *r, apop_model *m);
void apop_prep(apop_data *d, apop_model *m);
apop_model *apop_parameter_model(apop_data *d, apop_model *m);
apop_data * apop_predict(apop_data *d, apop_model *m);
apop_model *apop_beta_from_mean_var(double m, double v); //in apop_beta.c
#define apop_model_set_parameters(in, ...) apop_model_set_parameters_base((in), (double []) {__VA_ARGS__})
apop_model *apop_model_set_parameters_base(apop_model *in, double ap[]);
//apop_mixture.c
/** Produce a model as a linear combination of other models. See the documentation for the \ref apop_mixture model.
\param ... A list of models, either all parameterized or all unparameterized. See
examples in the \ref apop_mixture documentation.
*/
#define apop_model_mixture(...) apop_model_mixture_base((apop_model *[]){__VA_ARGS__, NULL})
apop_model *apop_model_mixture_base(apop_model **inlist);
//transform/apop_cross.c.
/** Generate a model consisting of the cross product of several independent models. The output \ref apop_model
is a copy of \ref apop_cross; see that model's documentation for details.
\li If you input only one model, return a copy of that model; print a warning iff <tt>apop_opts.verbose >= 2</tt>.
\exception error=='n' First model input is \c NULL.
Examples:
\include cross_models.c
*/
#define apop_model_cross(...) apop_model_cross_base((apop_model *[]){__VA_ARGS__, NULL})
apop_model *apop_model_cross_base(apop_model *mlist[]);
////More functions
//The variadic versions, with lots of options to input extra parameters to the
//function being mapped/applied
Apop_var_declare( apop_data * apop_map(apop_data *in, double (*fn_d)(double), double (*fn_v)(gsl_vector*),
double (*fn_r)(apop_data *), double (*fn_dp)(double! void *), double (*fn_vp)(gsl_vector*! void *),
double (*fn_rp)(apop_data *! void *), double (*fn_dpi)(double! void *! int),
double (*fn_vpi)(gsl_vector*! void *! int), double (*fn_rpi)(apop_data*! void *! int),
double (*fn_di)(double! int), double (*fn_vi)(gsl_vector*! int), double (*fn_ri)(apop_data*! int),
void *param, int inplace, char part, int all_pages) )
Apop_var_declare( double apop_map_sum(apop_data *in, double (*fn_d)(double), double (*fn_v)(gsl_vector*),
double (*fn_r)(apop_data *), double (*fn_dp)(double! void *), double (*fn_vp)(gsl_vector*! void *),
double (*fn_rp)(apop_data *! void *), double (*fn_dpi)(double! void *! int),
double (*fn_vpi)(gsl_vector*! void *! int), double (*fn_rpi)(apop_data*! void *! int),
double (*fn_di)(double! int), double (*fn_vi)(gsl_vector*! int), double (*fn_ri)(apop_data*! int),
void *param, char part, int all_pages) )
//the specific-to-a-type versions, quicker and easier when appropriate.
gsl_vector *apop_matrix_map(const gsl_matrix *m, double (*fn)(gsl_vector*));
gsl_vector *apop_vector_map(const gsl_vector *v, double (*fn)(double));
void apop_matrix_apply(gsl_matrix *m, void (*fn)(gsl_vector*));
void apop_vector_apply(gsl_vector *v, void (*fn)(double*));
gsl_matrix * apop_matrix_map_all(const gsl_matrix *in, double (*fn)(double));
void apop_matrix_apply_all(gsl_matrix *in, void (*fn)(double *));
double apop_vector_map_sum(const gsl_vector *in, double(*fn)(double));
double apop_matrix_map_sum(const gsl_matrix *in, double (*fn)(gsl_vector*));
double apop_matrix_map_all_sum(const gsl_matrix *in, double (*fn)(double));
// Some output routines
Apop_var_declare( void apop_matrix_print(const gsl_matrix *data, char const *output_name, FILE *output_pipe, char output_type, char output_append) )
Apop_var_declare( void apop_vector_print(gsl_vector *data, char const *output_name, FILE *output_pipe, char output_type, char output_append) )
Apop_var_declare( void apop_data_print(const apop_data *data, char const *output_name, FILE *output_pipe, char output_type, char output_append) )
void apop_matrix_show(const gsl_matrix *data);
void apop_vector_show(const gsl_vector *data);
void apop_data_show(const apop_data *data);
//statistics
Apop_var_declare( double apop_vector_mean(gsl_vector const *v, gsl_vector const *weights))
Apop_var_declare( double apop_vector_var(gsl_vector const *v, gsl_vector const *weights))
Apop_var_declare( double apop_vector_skew_pop(gsl_vector const *v, gsl_vector const *weights))
Apop_var_declare( double apop_vector_kurtosis_pop(gsl_vector const *v, gsl_vector const *weights))
Apop_var_declare( double apop_vector_cov(gsl_vector const *v1, gsl_vector const *v2,
gsl_vector const *weights))
Apop_var_declare( double apop_vector_distance(const gsl_vector *ina, const gsl_vector *inb, const char metric, const double norm) )
Apop_var_declare( void apop_vector_normalize(gsl_vector *in, gsl_vector **out, const char normalization_type) )
apop_data * apop_data_covariance(const apop_data *in);
apop_data * apop_data_correlation(const apop_data *in);
long double apop_vector_entropy(gsl_vector *in);
long double apop_matrix_sum(const gsl_matrix *m);
double apop_matrix_mean(const gsl_matrix *data);
void apop_matrix_mean_and_var(const gsl_matrix *data, double *mean, double *var);
apop_data * apop_data_summarize(apop_data *data);
Apop_var_declare( double * apop_vector_percentiles(gsl_vector *data, char rounding) )
apop_data *apop_test_fisher_exact(apop_data *intab); //in apop_fisher.c
//from apop_t_f_chi.c:
Apop_var_declare( int apop_matrix_is_positive_semidefinite(gsl_matrix *m, char semi) )
double apop_matrix_to_positive_semidefinite(gsl_matrix *m);
long double apop_multivariate_gamma(double a, int p);
long double apop_multivariate_lngamma(double a, int p);
//apop_tests.c
apop_data * apop_t_test(gsl_vector *a, gsl_vector *b);
apop_data * apop_paired_t_test(gsl_vector *a, gsl_vector *b);
Apop_var_declare( apop_data* apop_anova(char *table, char *data, char *grouping1, char *grouping2) )
#define apop_ANOVA apop_anova
Apop_var_declare( apop_data * apop_f_test (apop_model *est, apop_data *contrast) )
#define apop_F_test apop_f_test
//from the regression code:
#define apop_estimate_r_squared(in) apop_estimate_coefficient_of_determination(in)
apop_data * apop_text_unique_elements(const apop_data *d, size_t col);
gsl_vector * apop_vector_unique_elements(const gsl_vector *v);
Apop_var_declare( apop_data * apop_data_to_factors(apop_data *data, char intype, int incol, int outcol) )
Apop_var_declare( apop_data * apop_data_get_factor_names(apop_data *data, int col, char type) )
Apop_var_declare( apop_data * apop_data_to_dummies(apop_data *d, int col, char type, int keep_first, char append, char remove) )
Apop_var_declare( long double apop_model_entropy(apop_model *in, int draws) )
Apop_var_declare( long double apop_kl_divergence(apop_model *from, apop_model *to, int draw_ct, gsl_rng *rng) )
apop_data *apop_estimate_coefficient_of_determination (apop_model *);
void apop_estimate_parameter_tests (apop_model *est);
//Bootstrapping & RNG
apop_data * apop_jackknife_cov(apop_data *data, apop_model *model);
Apop_var_declare( apop_data * apop_bootstrap_cov(apop_data *data, apop_model *model, gsl_rng* rng, int iterations, char keep_boots, char ignore_nans, apop_data **boot_store) )
gsl_rng *apop_rng_alloc(int seed);
double apop_rng_GHgB3(gsl_rng * r, double* a); //in apop_asst.c
#define apop_rng_get_thread(thread_in) apop_rng_get_thread_base(#thread_in[0]=='\0' ? -1: (thread_in+0))
gsl_rng *apop_rng_get_thread_base(int thread);
int apop_arms_draw (double *out, gsl_rng *r, apop_model *m);
// maximum likelihod estimation related functions
Apop_var_declare( gsl_vector * apop_numerical_gradient(apop_data * data, apop_model* model, double delta) )
Apop_var_declare( apop_data * apop_model_hessian(apop_data * data, apop_model *model, double delta) )
Apop_var_declare( apop_data * apop_model_numerical_covariance(apop_data * data, apop_model *model, double delta) )
void apop_maximum_likelihood(apop_data * data, apop_model *dist);
Apop_var_declare( apop_model * apop_estimate_restart (apop_model *e, apop_model *copy, char * starting_pt, double boundary) )
//in apop_linear_constraint.c
Apop_var_declare( long double apop_linear_constraint(gsl_vector *beta, apop_data * constraint, double margin) )
//in apop_model_fix_params.c
apop_model * apop_model_fix_params(apop_model *model_in);
apop_model * apop_model_fix_params_get_base(apop_model *model_in);
//////vtables
/** \cond doxy_ignore */
/* This declares the vtable macros for each procedure that uses the mechanism.
--We want to have type-checking on the functions put into the vtables. Type checking
happens only with functions, not macros, so we need a type_check function for every
vtable.
--Only once in your codebase, you'll need to #define Declare_type_checking_fns to
actually define the type checking function. Everywhere else, the function is merely
declared.
--All other uses point to having a macro, such as using __VA_ARGS__ to allow any sort
of inputs to the hash.
--We want to have such a macro for every vtable. That means that we need a macro
to write macros. We can't do that with C macros, so this file uses m4 macros to
generate C macros.
--After the m4 definition of make_vtab_fns, each new vtable requires a typedef, a hash
definition, and a call to make_vtab_fns to do the rest.
*/
m4_define(make_vtab_fns, <|m4_dnl
#ifdef Declare_type_checking_fns
void $1_type_check($1_type in){ };
#else
void $1_type_check($1_type in);
#endif
#define $1_vtable_add(fn, ...) $1_type_check(fn), apop_vtable_add("$1", fn, $1_hash(__VA_ARGS__))
#define $1_vtable_get(...) apop_vtable_get("$1", $1_hash(__VA_ARGS__))
#define $1_vtable_drop(...) apop_vtable_drop("$1", $1_hash(__VA_ARGS__))m4_dnl
|>)
int apop_vtable_add(char const *tabname, void *fn_in, unsigned long hash);
void *apop_vtable_get(char const *tabname, unsigned long hash);
int apop_vtable_drop(char const *tabname, unsigned long hash);
typedef apop_model *(*apop_update_type)(apop_data *, apop_model* , apop_model*);
#define apop_update_hash(m1, m2) ( \
((m1)->log_likelihood ? (size_t)(m1)->log_likelihood : \
(m1)->p ? (size_t)(m1)->p*33 : \
(m1)->draw ? (size_t)(m1)->draw*33*27 \
: 33*27*19) \
+((m2)->log_likelihood ? (size_t)(m2)->log_likelihood : \
(m2)->p ? (size_t)(m2)->p*33 : \
(m2)->draw ? (size_t)(m2)->draw*33*27 \
: 33*27*19 \
) * 37)
make_vtab_fns(apop_update)
typedef long double (*apop_entropy_type)(apop_model *model);
#define apop_entropy_hash(m1) ((size_t)(m1)->log_likelihood + 33 * (size_t)((m1)->p) + 27*(size_t)((m1)->draw))
make_vtab_fns(apop_entropy)
typedef void (*apop_score_type)(apop_data *d, gsl_vector *gradient, apop_model *params);
#define apop_score_hash(m1) ((size_t)((m1)->log_likelihood ? (m1)->log_likelihood : (m1)->p))
make_vtab_fns(apop_score)
typedef apop_model* (*apop_parameter_model_type)(apop_data *, apop_model *);
#define apop_parameter_model_hash(m1) ((size_t)((m1)->log_likelihood ? (m1)->log_likelihood : (m1)->p)*33 + (m1)->estimate ? (size_t)(m1)->estimate: 27)
make_vtab_fns(apop_parameter_model)
typedef apop_data * (*apop_predict_type)(apop_data *d, apop_model *params);
#define apop_predict_hash(m1) ((size_t)((m1)->log_likelihood ? (m1)->log_likelihood : (m1)->p)*33 + (m1)->estimate ? (size_t)(m1)->estimate: 27)
make_vtab_fns(apop_predict)
typedef void (*apop_model_print_type)(apop_model *params, FILE *out);
#define apop_model_print_hash(m1) ((m1)->log_likelihood ? (size_t)(m1)->log_likelihood : \
(m1)->p ? (size_t)(m1)->p*33 : \
(m1)->estimate ? (size_t)(m1)->estimate*33*33 : \
(m1)->draw ? (size_t)(m1)->draw*33*27 : \
(m1)->cdf ? (size_t)(m1)->cdf*27*27 \
: 27)
make_vtab_fns(apop_model_print)
/** \endcond */ //End of Doxygen ignore.
//////Asst
long double apop_generalized_harmonic(int N, double s) __attribute__ ((__pure__));
apop_data * apop_test_anova_independence(apop_data *d);
#define apop_test_ANOVA_independence(d) apop_test_anova_independence(d)
Apop_var_declare( int apop_regex(const char *string, const char* regex, apop_data **substrings, const char use_case) )
int apop_system(const char *fmt, ...) __attribute__ ((format (printf,1,2)));
//Histograms and PMFs
gsl_vector * apop_vector_moving_average(gsl_vector *, size_t);
apop_data * apop_histograms_test_goodness_of_fit(apop_model *h0, apop_model *h1);
apop_data * apop_test_kolmogorov(apop_model *m1, apop_model *m2);
apop_data *apop_data_pmf_compress(apop_data *in);
Apop_var_declare( apop_data * apop_data_to_bins(apop_data const *indata, apop_data const *binspec, int bin_count, char close_top_bin) )
Apop_var_declare( apop_model * apop_model_to_pmf(apop_model *model, apop_data *binspec, long int draws, int bin_count) )
//text conveniences
Apop_var_declare( char* apop_text_paste(apop_data const*strings, char *between, char *before, char *after, char *between_cols, int (*prune)(apop_data* ! int ! int ! void*), void* prune_parameter) )
/** Notify the user of errors, warning, or debug info.
writes to \ref apop_opts.log_file, which is a \c FILE handle. The default is \c stderr,
but use \c fopen to attach to a file.
\param verbosity At what verbosity level should the user be warned? E.g., if level==2, then print iff apop_opts.verbosity >= 2.
\param ... The message to write to the log (presuming the verbosity level is high
enough). This can be a printf-style format with following arguments,
e.g., <tt>apop_notify(0, "Beta is currently %g", beta)</tt>.
*/
#define Apop_notify(verbosity, ...) {\
if (apop_opts.verbose != -1 && apop_opts.verbose >= verbosity) { \
if (!apop_opts.log_file) apop_opts.log_file = stderr; \
fprintf(apop_opts.log_file, "%s: ", __func__); fprintf(apop_opts.log_file, __VA_ARGS__); fprintf(apop_opts.log_file, "\n"); \
fflush(apop_opts.log_file); \
} }
/** \cond doxy_ignore */
#define Apop_maybe_abort(level) \
{if ((apop_opts.verbose >= level && apop_opts.stop_on_warning == 'v') \
|| (apop_opts.stop_on_warning=='w') ) \
raise(SIGTRAP);}
/** \endcond */
/** Execute an action and print a message to the current \c FILE handle held by <tt>apop_opts.log_file</tt> (default: \c stderr).
\param test The expression that, if true, triggers the action.
\param onfail If the assertion fails, do this. E.g., <tt>out->error='x'; return GSL_NAN</tt>. Notice that it is OK to include several lines of semicolon-separated code here, but if you have a lot to do, the most readable option may be <tt>goto outro</tt>, plus an appropriately-labeled section at the end of your function.
\param level Print the warning message only if \ref apop_opts_type "apop_opts.verbose" is greater than or equal to this. Zero usually works, but for minor infractions use one, or for more verbose debugging output use 2.
\param ... The error message in printf form, plus any arguments to be inserted into the printf string. I'll provide the function name and a carriage return.
Some examples:
\code
//the typical case, stopping function execution:
Apop_stopif(isnan(x), return NAN, 0, "x is NAN; failing");
//Mark a flag, go to a cleanup step
Apop_stopif(x < 0, needs_cleanup=1; goto cleanup, 0, "x is %g; cleaning up and exiting.", x);
//Print a diagnostic iff <tt>apop_opts.verbose>=1</tt> and continue
Apop_stopif(x < 0, , 1, "warning: x is %g.", x);
\endcode
\li If \c apop_opts.stop_on_warning is nonzero and not <tt>'v'</tt>, then a failed test halts via \c abort(), even if the <tt>apop_opts.verbose</tt> level is set so that the warning message doesn't print to screen. Use this when running via debugger.
\li If \c apop_opts.stop_on_warning is <tt>'v'</tt>, then a failed test halts via \c abort() iff the verbosity level is high enough to print the error.
*/
#define Apop_stopif(test, onfail, level, ...) do {\
if (test) { \
Apop_notify(level, __VA_ARGS__); \
Apop_maybe_abort(level) \
onfail; \
} } while(0)
#define apop_errorlevel -5
/** \cond doxy_ignore */
//For use in stopif, to return a blank apop_data set with an error attached.
#define apop_return_data_error(E) {apop_data *out=apop_data_alloc(); out->error='E'; return out;}
/* The Apop_stopif macro is currently favored, but there's a long history of prior
error-handling setups. Consider all of the Assert... macros below to be deprecated.
*/
#define Apop_assert_c(test, returnval, level, ...) \
Apop_stopif(!(test), return returnval, level, __VA_ARGS__)
#define Apop_assert(test, ...) Apop_assert_c((test), 0, apop_errorlevel, __VA_ARGS__)
//For things that return void. Transitional and deprecated at birth.
#define Apop_assert_n(test, ...) Apop_assert_c((test), , apop_errorlevel, __VA_ARGS__)
#define Apop_assert_negone(test, ...) Apop_assert_c((test), -1, apop_errorlevel, __VA_ARGS__)
/** \endcond */ //End of Doxygen ignore.
//Missing data
Apop_var_declare( apop_data * apop_data_listwise_delete(apop_data *d, char inplace) )
apop_model * apop_ml_impute(apop_data *d, apop_model* meanvar);
Apop_var_declare(apop_model *apop_model_metropolis(apop_data *d, gsl_rng* rng, apop_model *m))
Apop_var_declare( apop_model * apop_update(apop_data *data, apop_model *prior, apop_model *likelihood, gsl_rng *rng) )
Apop_var_declare( double apop_test(double statistic, char *distribution, double p1, double p2, char tail) )
//apop_sort.c
Apop_var_declare( apop_data *apop_data_sort(apop_data *data, apop_data *sort_order, char asc, char inplace, double *col_order))
//raking
Apop_var_declare( apop_data * apop_rake(char const *margin_table, char * const*var_list,
int var_ct, char * const *contrasts, int contrast_ct,
char const *structural_zeros, int max_iterations, double tolerance,
char const *count_col, char const *init_table,
char const *init_count_col, double nudge) )
#include <gsl/gsl_cdf.h>
#include <gsl/gsl_blas.h>
#include <gsl/gsl_sf_log.h>
#include <gsl/gsl_sf_exp.h>
#include <gsl/gsl_linalg.h>
#include <gsl/gsl_sf_gamma.h>
#include <gsl/gsl_sf_psi.h>
#include <gsl/gsl_randist.h>
#include <gsl/gsl_histogram.h>
#include <gsl/gsl_statistics_double.h>
//Some linear algebra utilities
double apop_det_and_inv(const gsl_matrix *in, gsl_matrix **out, int calc_det, int calc_inv);
Apop_var_declare( apop_data * apop_dot(const apop_data *d1, const apop_data *d2, char form1, char form2) )
Apop_var_declare( int apop_vector_bounded(const gsl_vector *in, long double max) )
gsl_matrix * apop_matrix_inverse(const gsl_matrix *in) ;
double apop_matrix_determinant(const gsl_matrix *in) ;
//apop_data* apop_sv_decomposition(gsl_matrix *data, int dimensions_we_want);
Apop_var_declare( apop_data * apop_matrix_pca(gsl_matrix *data, int const dimensions_we_want) )
Apop_var_declare( gsl_vector * apop_vector_stack(gsl_vector *v1, gsl_vector const * v2, char inplace) )
Apop_var_declare( gsl_matrix * apop_matrix_stack(gsl_matrix *m1, gsl_matrix const * m2, char posn, char inplace) )
void apop_vector_log(gsl_vector *v);
void apop_vector_log10(gsl_vector *v);
void apop_vector_exp(gsl_vector *v);
////Subsetting macros
/** \cond doxy_ignore */
/** These are all deprecated.*/
#define APOP_SUBMATRIX(m, srow, scol, nrows, ncols, o) gsl_matrix apop_mm_##o = gsl_matrix_submatrix((m), (srow), (scol), (nrows),(ncols)).matrix;\
gsl_matrix * o = &( apop_mm_##o ); // Use \ref Apop_subm.
#define Apop_submatrix APOP_SUBMATRIX
#define Apop_col_v(m, col, v) gsl_vector apop_vv_##v = ((col) == -1) ? (gsl_vector){} : gsl_matrix_column((m)->matrix, (col)).vector;\
gsl_vector * v = ((col)==-1) ? (m)->vector : &( apop_vv_##v ); // Use \ref Apop_cv.
#define Apop_row_v(m, row, v) Apop_matrix_row((m)->matrix, row, v) // Use \ref Apop_rv.
#define Apop_rows(d, rownum, len, outd) apop_data *outd = Apop_rs(d, rownum, len) // Use \ref Apop_rs.
#define Apop_row(d, row, outd) Apop_rows(d, row, 1, outd) // Use \ref Apop_r.
#define Apop_cols(d, colnum, len, outd) apop_data *outd = Apop_cs(d, colnum, len); // Use \ref Apop_cs.
/** \endcond */ //End of Doxygen ignore.
/** \def Apop_row_tv(m, row_name, v)
After this call, \c v will hold a \c gsl_vector view of an \ref apop_data set \c m. The view will consist only of the row with name \c row_name.
Unlike \ref Apop_rv, the second argument is a row name, that I'll look up using \ref apop_name_find, and the third is the name of the view to be generated.
\see Apop_rs, Apop_r, Apop_rv, Apop_row_t, Apop_mrv
*/
#define Apop_row_tv(m, row, v) gsl_vector apop_vv_##v = gsl_matrix_row((m)->matrix, apop_name_find((m)->names, row, 'r')).vector;\
gsl_vector * v = &( apop_vv_##v );
/** \def Apop_col_tv(m, col_name, v)
After this call, \c v will hold a \c gsl_vector view of the \ref apop_data set \c m.
The view will consist only of the column with name \c col_name.
Unlike \ref Apop_cv, the second argument is a column name, that I'll look up using \ref apop_name_find, and the third is the name of the view to be generated.
\see Apop_cs, Apop_c, Apop_cv, Apop_col_t, Apop_mcv
*/
#define Apop_col_tv(m, col, v) gsl_vector apop_vv_##v = gsl_matrix_column((m)->matrix, apop_name_find((m)->names, col, 'c')).vector;\
gsl_vector * v = &( apop_vv_##v );
/** \def Apop_row_t(m, row_name, v)
After this call, \c v will hold an \ref apop_data view of an \ref apop_data set \c m. The view will consist only of the row with name \c row_name.
Unlike \ref Apop_r, the second argument is a row name, that I'll look up using \ref apop_name_find, and the third is the name of the view to be generated.
\see Apop_rs, Apop_r, Apop_rv, Apop_row_tv, Apop_mrv
*/
#define Apop_row_t(d, rowname, outd) int apop_row_##outd = apop_name_find((d)->names, rowname, 'r'); Apop_rows(d, apop_row_##outd, 1, outd)
/** \def Apop_col_t(m, col_name, v)
After this call, \c v will hold a view of the \ref apop_data set \c m. The view will consist only of a \c gsl_vector view of the column of the \ref apop_data set \c m with name \c col_name.
Unlike \ref Apop_c, the second argument is a column name, that I'll look up using \ref apop_name_find, and the third is the name of the view to be generated.
\see Apop_cs, Apop_c, Apop_cv, Apop_col_tv, Apop_mcv
*/
#define Apop_col_t(d, colname, outd) int apop_col_##outd = apop_name_find((d)->names, colname, 'c'); Apop_cols(d, apop_col_##outd, 1, outd)
// The above versions relied on gsl_views, which stick to C as of 1989 CE.
// Better to just create the views via designated initializers.
/** \def Apop_subm(data_to_view, srow, scol, nrows, ncols)
Generate a view of a submatrix within a \c gsl_matrix. Like \ref Apop_r, et al., the view is an automatically-allocated variable that is lost once the program flow leaves the scope in which it is declared.
\param data_to_view The root matrix
\param srow the first row (in the root matrix) of the top of the submatrix
\param scol the first column (in the root matrix) of the left edge of the submatrix
\param nrows number of rows in the submatrix
\param ncols number of columns in the submatrix
\return An automatically-allocated view of type \c gsl_matrix.
*/
#define Apop_subm(matrix_to_view, srow, scol, nrows, ncols)( \
(!(matrix_to_view) \
|| (matrix_to_view)->size1 < (srow)+(nrows) || (srow) < 0 \
|| (matrix_to_view)->size2 < (scol)+(ncols) || (scol) < 0) ? NULL \
: &(gsl_matrix){.size1=(nrows), .size2=(ncols), \
.tda=(matrix_to_view)->tda, \
.data=gsl_matrix_ptr((matrix_to_view), (srow), (scol))} \
)
/** Get a vector view of a single row of a \ref gsl_matrix.
\param matrix_to_vew A \ref gsl_matrix.
\param row An integer giving the row to be viewed.
\return A \c gsl_vector view of the given row. The view is automatically allocated,
and disappears as soon as the program leaves the scope in which it is declared.
See \ref apop_vector_correlation for an example of use.
\see Apop_r, Apop_rv
*/
#define Apop_mrv(matrix_to_view, row) Apop_rv(&(apop_data){.matrix=matrix_to_view}, row)
/** Get a vector view of a single column of a \ref gsl_matrix.
\param matrix_to_vew A \ref gsl_matrix.
\param row An integer giving the column to be viewed.
\return A \c gsl_vector view of the given column. The view is automatically allocated,
and disappears as soon as the program leaves the scope in which it is declared.
\code
gsl_matrix *m = apop_query_to_data("select col1, col2, col3 from data")->matrix;
printf("The correlation coefficient between columns two "
"and three is %g.\n", apop_vector_correlation(Apop_mcv(m, 2), Apop_mcv(m, 3)));
\endcode
\see Apop_r, Apop_cv
*/
#define Apop_mcv(matrix_to_view, col) Apop_cv(&(apop_data){.matrix=matrix_to_view}, col)
/** \def Apop_rv(d, row)
A macro to generate a temporary one-row view of the matrix in an \ref apop_data set \c d, pulling out only
row \c row. The view is a \c gsl_vector set.
\code
gsl_vector *v = Apop_rv(your_data, i);
for (int i=0; i< your_data->matrix->size1; i++)
printf("Σ_%i = %g\n", i, apop_vector_sum(Apop_r(your_data, i)));
\endcode
The view is automatically allocated, and disappears as soon as the program leaves the scope in which it is declared.
\see Apop_r, Apop_rv, Apop_row_tv, Apop_row_t, Apop_mrv
*/
#define Apop_rv(data_to_view, row) ( \
((data_to_view) == NULL || (data_to_view)->matrix == NULL \
|| (data_to_view)->matrix->size1 <= (row) || (row) < 0) ? NULL \
: &(gsl_vector){.size=(data_to_view)->matrix->size2, \
.stride=1, .data=gsl_matrix_ptr((data_to_view)->matrix, (row), 0)} \
)
/** \def Apop_cv(d, col)
A macro to generate a temporary one-column view of the matrix in an \ref apop_data
set \c d, pulling out only column \c col. The view is a \c gsl_vector set.
As usual, column -1 is the vector element of the \ref apop_data set.
\code
gsl_vector *v = Apop_cv(your_data, i);
for (int i=0; i< your_data->matrix->size2; i++)
printf("Σ_%i = %g\n", i, apop_vector_sum(Apop_c(your_data, i)));
\endcode
The view is automatically allocated, and disappears as soon as the program leaves the
scope in which it is declared.
\see Apop_cs, Apop_c, Apop_col_tv, Apop_col_t, Apop_mcv
*/
#define Apop_cv(data_to_view, col) ( \
!(data_to_view) ? NULL \
: (col)==-1 ? (data_to_view)->vector \
: (!(data_to_view)->matrix \
|| (data_to_view)->matrix->size2 <= (col) || ((int)(col)) < -1) ? NULL \
: &(gsl_vector){.size=(data_to_view)->matrix->size1, \
.stride=(data_to_view)->matrix->tda, .data=gsl_matrix_ptr((data_to_view)->matrix, 0, (col))} \
)
/** \cond doxy_ignore */
/* Not (yet) for public use. */
#define Apop_subvector(v, start, len) ( \
((v) == NULL || (v)->size < ((start)+(len)) || (start) < 0) ? NULL \
: &(gsl_vector){.size=(len), .stride=(v)->stride, .data=(v)->data+(start*(v)->stride)})
/** \endcond */
/** \def Apop_rs(d, row, len)
A macro to generate a temporary view of \ref apop_data set \c d pulling only certain rows, beginning at row \c row
and having height \c len.
The view is automatically allocated, and disappears as soon as the program leaves the scope in which it is declared.
\see Apop_r, Apop_rv, Apop_row_tv, Apop_row_t, Apop_mrv
*/
#define Apop_rs(d, rownum, len)( \
(!(d) || (rownum) < 0) ? NULL \
: &(apop_data){ \
.names= ( !((d)->names) ? NULL : \
&(apop_name){ \
.title = (d)->names->title, \
.vector = (d)->names->vector, \
.col = (d)->names->col, \
.row = ((d)->names->row && (d)->names->rowct > (rownum)) ? &((d)->names->row[rownum]) : NULL, \
.texthash = (d)->names->texthash, \
.rowhash = ((d)->names->rowhash && (d)->names->rowct > (rownum)) ? &((d)->names->rowhash[rownum]) : NULL, \
.colhash = (d)->names->colhash, \
.text = (d)->names->text, \
.colct = (d)->names->colct, \
.rowct = (d)->names->row ? (GSL_MIN(1, GSL_MAX((d)->names->rowct - (int)(rownum), 0))) \
: 0, \
.textct = (d)->names->textct }), \
.vector= Apop_subvector((d->vector), (rownum), (len)), \
.matrix = Apop_subm(((d)->matrix), (rownum), 0, (len), (d)->matrix?(d)->matrix->size2:0), \
.weights = Apop_subvector(((d)->weights), (rownum), (len)), \
.textsize[0]=(d)->textsize[0]> (rownum)+(len)-1 ? (len) : 0, \
.textsize[1]=(d)->textsize[1], \
.text = (d)->text ? &((d)->text[rownum]) : NULL, \
})
/** \def Apop_cs(d, col, len)
A macro to generate a temporary view of \ref apop_data set \c d including only certain columns, beginning at column \c col and having length \c len.
The view is automatically allocated, and disappears as soon as the program leaves the scope in which it is declared.
\see Apop_c, Apop_cv, Apop_col_tv, Apop_col_t, Apop_mcv
*/
#define Apop_cs(d, colnum, len) ( \
(!(d)||!(d)->matrix || (d)->matrix->size2 <= (colnum)+(len)-1) \
? NULL \
: &(apop_data){ \
.vector= NULL, \
.weights= (d)->weights, \
.matrix = Apop_subm((d)->matrix, 0, colnum, (d)->matrix->size1, (len)),\
.textsize[0] = 0, \
.textsize[1] = 0, \
.text = NULL, \
.names= (d)->names ? &(apop_name){ \
.title = (d)->names->title, \
.vector = NULL, \
.row = (d)->names->row, \
.col = ((d)->names->col && (d)->names->colct > colnum) ? &((d)->names->col[colnum]) : NULL, \
.text = NULL, \
.texthash = NULL, \
.rowhash = (d)->names->rowhash, \
.colhash = ((d)->names->colhash && (d)->names->colct > (colnum)) ? &((d)->names->colhash[colnum]) : NULL, \
.rowct = (d)->names->rowct, \
.colct = (d)->names->col ? (GSL_MIN(len, GSL_MAX((d)->names->colct - colnum, 0))) \
: 0, \
.textct = (d)->names->textct } : NULL \
})
/** \def Apop_r(d, row)
A macro to generate a temporary one-row view of \ref apop_data set \c d, pulling out only
row \c row. The view is also an \ref apop_data set, with names and other decorations.
\code
//pull a single row
apop_data *v = Apop_r(your_data, 7);
//or loop through a sequence of one-row data sets.
apop_model *std = apop_model_set_parameters(apop_normal, 0, 1);
for (int i=0; i< your_data->matrix->size1; i++)
printf("Std Normal CDF up to observation %i is %g\n",
i, apop_cdf(Apop_r(your_data, i), std));
\endcode
The view is automatically allocated, and disappears as soon as the program leaves the
scope in which it is declared.
\see Apop_rs, Apop_row_v, Apop_row_tv, Apop_row_t, Apop_mrv
*/
#define Apop_r(d, rownum) Apop_rs(d, rownum, 1)
/** \def Apop_c(d, col)
A macro to generate a temporary one-column view of \ref apop_data set \c d, pulling out only
column \c col.
After this call, \c outd will be a pointer to this temporary
view, that you can use as you would any \ref apop_data set.
\see Apop_cs, Apop_cv, Apop_col_tv, Apop_col_t, Apop_mcv
*/
#define Apop_c(d, col) Apop_cs(d, col, 1)
/** \cond doxy_ignore */
#define APOP_COL Apop_col
#define apop_col Apop_col
#define APOP_COL_T Apop_col_t
#define apop_col_t Apop_col_t
#define APOP_COL_TV Apop_col_tv
#define apop_col_tv Apop_col_tv
#define APOP_ROW Apop_row
#define apop_row Apop_row
#define APOP_COLS Apop_cols
#define apop_cols Apop_cols
#define APOP_COL_V Apop_col_v
#define apop_col_v Apop_col_v
#define APOP_ROW_V Apop_row_v
#define apop_row_v Apop_row_v
#define APOP_ROWS Apop_rows
#define apop_rows Apop_rows
#define Apop_data_row Apop_row #deprecated
#define APOP_ROW_T Apop_row_t
#define apop_row_t Apop_row_t
#define APOP_ROW_TV Apop_row_tv
#define apop_row_tv Apop_row_tv
/** Deprecated. Use Apop_mrv */
#define Apop_matrix_row(m, row, v) gsl_vector apop_vv_##v = gsl_matrix_row((m), (row)).vector;\
gsl_vector * v = &( apop_vv_##v );
/* Deprecated. Use Apop_mcv */
#define Apop_matrix_col(m, col, v) gsl_vector apop_vv_##v = gsl_matrix_column((m), (col)).vector;\
gsl_vector * v = &( apop_vv_##v );
#define APOP_MATRIX_ROW Apop_matrix_row
#define apop_matrix_row Apop_matrix_row
#define APOP_MATRIX_COL Apop_matrix_col
#define apop_matrix_col Apop_matrix_col
/** \endcond */
long double apop_vector_sum(const gsl_vector *in);
double apop_vector_var_m(const gsl_vector *in, const double mean);
Apop_var_declare( double apop_vector_correlation(const gsl_vector *ina, const gsl_vector *inb, const gsl_vector *weights) )
double apop_vector_kurtosis(const gsl_vector *in);
double apop_vector_skew(const gsl_vector *in);
#define apop_sum apop_vector_sum
#define apop_var apop_vector_var
#define apop_mean apop_vector_mean
//////database utilities
Apop_var_declare( int apop_table_exists(char const *name, char remove) )
int apop_db_open(char const *filename);
Apop_var_declare( int apop_db_close(char vacuum) )