forked from b-k/apophenia
-
Notifications
You must be signed in to change notification settings - Fork 0
/
apop_data.m4.c
1530 lines (1374 loc) · 74.2 KB
/
apop_data.m4.c
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
The apop_data structure joins together a gsl_matrix, apop_name, and a table of strings. */
/* Copyright (c) 2006--2009 by Ben Klemens. Licensed under the modified GNU GPL v2; see COPYING and COPYING2. */
#include "apop_internal.h"
//apop_gsl_error is in apop_linear_algebra.c
#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);
/** Allocate a \ref apop_data structure, to be filled with data.
\li The typical case is three arguments, like <tt>apop_data_alloc(2,3,4)</tt>: vector size, matrix rows, matrix cols. If the first argument is zero, you get a \c NULL vector.
\li Two arguments, <tt>apop_data_alloc(2,3)</tt>, would allocate just a matrix, leaving the vector \c NULL.
\li One argument, <tt>apop_data_alloc(2)</tt>, would allocate just a vector, leaving the matrix \c NULL.
\li Zero arguments, <tt>apop_data_alloc()</tt>, will produce a basically blank set, with \c out->matrix==out->vector==NULL.
For allocating the text part, see \ref apop_text_alloc.
The \c weights vector is set to \c NULL. If you need it, allocate it via
\code d->weights = gsl_vector_alloc(row_ct); \endcode
\see apop_data_calloc
\return The \ref apop_data structure, allocated and ready.
\exception out->error=='a' Allocation error. The matrix, vector, or names couldn't be <tt>malloc</tt>ed, which probably means that you requested a very large data set.
\li An \ref apop_data struct, by itself, is about 72 bytes. If I can't allocate that much memory, I return \c NULL.
But if even this much fails, your computer may be on fire and you should go put it out.
\li This function uses the \ref designated syntax for inputs.
\ingroup data_struct
*/
APOP_VAR_HEAD apop_data * apop_data_alloc(const size_t size1, const size_t size2, const int size3){
const size_t apop_varad_var(size1, 0);
const size_t apop_varad_var(size2, 0);
const int apop_varad_var(size3, 0);
APOP_VAR_ENDHEAD
size_t vsize=0, msize1=0;
int msize2=0;
if (size3){
vsize = size1;
msize1 = size2;
msize2 = size3;
}
else if (size2) {
msize1 = size1;
msize2 = size2;
}
else vsize = size1;
apop_data *setme = malloc(sizeof(apop_data));
Apop_stopif(!setme, return NULL, -5, "malloc failed. Probably out of memory.");
*setme = (apop_data) { }; //init to zero/NULL.
Set_gsl_handler
if (msize2 > 0 && msize1 > 0){
setme->matrix = gsl_matrix_alloc(msize1,msize2);
Apop_stopif(!setme->matrix, setme->error='a'; return setme,
0, "malloc failed on a %zu x %i matrix. Probably out of memory.", msize1, msize2);
}
if (vsize){
setme->vector = gsl_vector_alloc(vsize);
Apop_stopif(!setme->vector, setme->error='a'; return setme,
0, "malloc failed on a vector of size %zu. Probably out of memory.", vsize);
}
Unset_gsl_handler
setme->names = apop_name_alloc();
Apop_stopif(!setme->names, setme->error='a'; return setme,
0, "couldn't allocate names. Probably out of memory.");
return setme;
}
/** Allocate a \ref apop_data structure, to be filled with data; set everything in the allocated portion to zero. See \ref apop_data_alloc for details.
\return The \ref apop_data structure, allocated and zeroed out.
\exception out->error=='m' malloc error; probably out of memory.
\see apop_data_alloc
\ingroup data_struct
\li This function uses the \ref designated syntax for inputs.
*/
APOP_VAR_HEAD apop_data * apop_data_calloc(const size_t size1, const size_t size2, const int size3){
const size_t apop_varad_var(size1, 0);
const size_t apop_varad_var(size2, 0);
const int apop_varad_var(size3, 0);
APOP_VAR_ENDHEAD
size_t vsize=0, msize1=0;
int msize2=0;
if (size3){
vsize = size1;
msize1 = size2;
msize2 = size3;
}
else if (size2) {
msize1 = size1;
msize2 = size2;
}
else vsize = size1;
apop_data *setme = malloc(sizeof(apop_data));
Apop_stopif(!setme, apop_return_data_error('m'), 0, "malloc failed. Probably out of memory.");
*setme = (apop_data) { }; //init to zero/NULL.
if (msize2 >0 && msize1 > 0){
setme->matrix = gsl_matrix_calloc(msize1,msize2);
Apop_stopif(!setme->matrix, apop_return_data_error('m'), 0, "malloc failed on a %zu x %i matrix. Probably out of memory.", msize1, msize2);
}
if (vsize){
setme->vector = gsl_vector_calloc(vsize);
Apop_stopif(!setme->vector, apop_return_data_error('m'), 0, "malloc failed on a vector of size %zu. Probably out of memory.", vsize);
}
setme->names = apop_name_alloc();
return setme;
}
/*For a touch of space saving, blank strings in a text grid
all point to the same nul string. */
char *apop_nul_string = "";
static void apop_text_blank(apop_data *in, const size_t row, const size_t col){
if (in->text[row][col] != apop_nul_string) free(in->text[row][col]);
in->text[row][col] = apop_nul_string;
}
/** Free a matrix of chars* (i.e., a char***). This is the form of the
text element of the \ref apop_data set, so you can use this for:
\code
apop_text_free(yourdata->text, yourdata->textsize[0], yourdata->textsize[1]);
\endcode
This is what \c apop_data_free uses internally.
*/
void apop_text_free(char ***freeme, int rows, int cols){
if (rows && cols)
for (int i=0; i < rows; i++){
for (int j=0; j < cols; j++)
if(freeme[i][j]!=apop_nul_string)
free(freeme[i][j]);
free(freeme[i]);
}
free(freeme);
}
/** Free the elements of the given \ref apop_data set and then the \ref apop_data set
itself. Intended to be used by \ref apop_data_free, a macro that calls this to free
elements, then sets the value to \c NULL.
\li \ref apop_data_free is a macro that calls this function and, on success, sets the input pointer to \c NULL.
For typical cases, that's slightly more useful than this function.
\exception freeme.error='c' Circular linking is against the rules. If <tt>freeme->more == freeme</tt>, then
I set <tt>freeme.error='c'</tt> and return. If you send in a structure like A -> B ->
B, then both data sets A and B will be marked.
\return \c 0 on OK, \c 'c' on error.
*/
char apop_data_free_base(apop_data *freeme){
if (!freeme) return 0;
if (freeme->more){
Apop_stopif(freeme == freeme->more, freeme->error='c'; return 'c',
1, "the ->more element of this data set equals the data set itself. "
"This is not healthy. Not freeing; marking your data set with error='c'.");
if (apop_data_free_base(freeme->more))
Apop_stopif(freeme->more->error == 'c', freeme->error='c'; return 'c',
1, "Propogating error code to parent data set");
}
if (freeme->vector)
gsl_vector_free(freeme->vector);
if (freeme->matrix)
gsl_matrix_free(freeme->matrix);
if (freeme->weights)
gsl_vector_free(freeme->weights);
apop_name_free(freeme->names);
apop_text_free(freeme->text, freeme->textsize[0] , freeme->textsize[1]);
free(freeme);
return 0;
}
/** Copy one \ref 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.
\li If you want space allocated or are unsure about dimensions, use \ref apop_data_copy.
\li If both \c in and \c out have a \c more pointer, also copy subsequent page(s).
\li You can use the subsetting macros, \ref Apop_r or \ref Apop_rs, to copy within a data set:
\code
//Copy the contents of row i of mydata to row j.
apop_data *fromrow = Apop_r(mydata, i);
apop_data *torow = Apop_r(mydata, j);
apop_data_memcpy(torow, fromrow);
// or just
apop_data_memcpy(Apop_r(mydata, i), Apop_r(mydata, j));
\endcode
\param out A structure that this function will fill. Must be preallocated with the appropriate sizes.
\param in The input data.
\exception out.error='d' Dimension error; couldn't copy.
\exception out.error='p' Part missing; e.g., in->matrix exists but out->matrix doesn't; couldn't copy.
\ingroup data_struct
*/
void apop_data_memcpy(apop_data *out, const apop_data *in){
Apop_stopif(!out, return, 0, "you are copying to a NULL matrix. Do you mean to use apop_data_copy instead?");
Apop_stopif(out==in, return, 1, "out==in. Doing nothing.");
if (in->matrix){
Apop_stopif(!out->matrix, out->error='p'; return, 1, "in->matrix exists but out->matrix does not.");
Apop_stopif(in->matrix->size1 != out->matrix->size1 || in->matrix->size2 != out->matrix->size2,
out->error='d'; return,
1, "you're trying to copy a (%zu X %zu) into a (%zu X %zu) matrix.",
in->matrix->size1, in->matrix->size2, out->matrix->size1, out->matrix->size2);
gsl_matrix_memcpy(out->matrix, in->matrix);
}
if (in->vector){
Apop_stopif(!out->vector, out->error='p'; return, 1, "in->vector exists but out->vector does not.");
Apop_stopif(in->vector->size != out->vector->size,
out->error='d'; return,
1, "You're trying to copy a %zu-elmt "
"vector into a %zu-elmt vector.", in->vector->size, out->vector->size);
gsl_vector_memcpy(out->vector, in->vector);
}
if (in->weights){
Apop_stopif(!out->weights, out->error='p'; return, 1, "in->weights exists but out->weights does not.");
Apop_stopif(in->weights->size != out->weights->size,
out->error='d'; return,
1, "Weight vector sizes don't match: "
"you're trying to copy a %zu-elmt vector into a %zu-elmt vector.",
in->weights->size, out->weights->size);
gsl_vector_memcpy(out->weights, in->weights);
}
if (in->names){
apop_name_free(out->names);
out->names = apop_name_alloc();
apop_name_stack(out->names, in->names, 'v');
apop_name_stack(out->names, in->names, 'r');
apop_name_stack(out->names, in->names, 'c');
apop_name_stack(out->names, in->names, 't');
}
out->textsize[0] = in->textsize[0];
out->textsize[1] = in->textsize[1];
if (in->textsize[0] && in->textsize[1]){
Apop_stopif(out->textsize[0] < in->textsize[0] || out->textsize[1] < in->textsize[1],
out->error='d'; return,
1, "I am trying to copy a grid of (%zu, %zu) text elements into a grid of (%zu, %zu), "
"and that won't work. Please use apop_text_alloc to reallocate the right amount of data, "
"or use apop_data_copy for automatic allocation.",
in->textsize[0] , in->textsize[1] , out->textsize[0] , out->textsize[1]);
for (size_t i=0; i< in->textsize[0]; i++)
for(size_t j=0; j < in->textsize[1]; j ++)
if (in->text[i][j] == apop_nul_string)
apop_text_blank(out, i, j);
else apop_text_add(out, i, j, "%s", in->text[i][j]);
}
if (in->more && out->more) apop_data_memcpy(out->more, in->more);
}
/** Copy one \ref apop_data structure to another. That is, all data is duplicated.
Basically a front-end for \ref apop_data_memcpy for those who prefer this sort of syntax.
Unlike \ref apop_data_memcpy, I do follow the \c more pointer.
\param in the input data
\return a structure that this function will allocate and fill. If input is NULL, then this will be NULL.
\ingroup data_struct
\exception out.error='a' Allocation error.
\exception out.error='c' Cyclic link: <tt>D->more == D</tt> (may be later in the chain, e.g., <tt>D->more->more = D->more</tt>) You'll have only a partial copy.
\exception out.error='d' Dimension error; should never happen.
\exception out.error='p' Missing part error; should never happen.
\li 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 <tt>>=1</tt>.
*/
apop_data *apop_data_copy(const apop_data *in){
if (!in) return NULL;
apop_data *out = apop_data_alloc();
Apop_stopif(out->error, return out, 0, "Allocation error.");
if (in->error){
Apop_notify(1, "the data set to be copied has an error flag of %c. Copying it.", in->error);
out->error = in->error;
}
if (in->more){
Apop_stopif(in == in->more, out->error='c'; return out,
0, "the ->more element of this data set equals the "
"data set itself. This is not healthy. Made a partial copy and set out.error='c'.");
out->more = apop_data_copy(in->more);
Apop_stopif(out->more->error, out->error=out->more->error; return out,
0, "propagating an error in the ->more element to the parent apop_data set. Only a partial copy made.");
}
if (in->vector){
out->vector = gsl_vector_alloc(in->vector->size);
Apop_stopif(!out->vector, out->error='a'; return out, 0, "Allocation error on vector of size %zu.", in->vector->size);
}
if (in->matrix){
out->matrix = gsl_matrix_alloc(in->matrix->size1, in->matrix->size2);
Apop_stopif(!out->matrix, out->error='a'; return out, 0, "Allocation error on matrix "
"of size %zu X %zu.", in->matrix->size1, in->matrix->size2);
}
if (in->weights){
out->weights = gsl_vector_alloc(in->weights->size);
Apop_stopif(!out->weights, out->error='a'; return out, 0, "Allocation error on weights vector of size %zu.", in->weights->size);
}
if (in->textsize[0] && in->textsize[1]){
apop_text_alloc(out, in->textsize[0], in->textsize[1]);
Apop_stopif(out->error, return out, 0, "Allocation error on text grid of size %zu X %zu.", in->textsize[0], in->textsize[1]);
}
apop_data_memcpy(out, in);
return out;
}
/** 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 \ref apop_data_split.
\param m1 the upper/rightmost data set (default = \c NULL)
\param m2 the second data set (default = \c NULL)
\param posn If 'r', stack rows of m1's matrix above rows of m2's<br>
if 'c', stack columns of m1's matrix to left of m2's<br>
(default = 'r')
\param inplace If \c 'y', use \ref apop_matrix_realloc and \ref apop_vector_realloc to modify \c m1 in place; see the caveats on those function. Otherwise, allocate a new vector, leaving \c m1 unmolested. (default='n')
\return The stacked data, either in a new \ref apop_data set or \c m1
\exception out->error=='a' Allocation error.
\exception out->error=='d' Dimension error; couldn't make a complete copy.
\li 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 \c m2 is \c NULL and \c inplace is \c 'y', where you'll get the original \c m1 pointer back)
\li 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].
\li \c more is ignored.
\li 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.
\li The same rules for dealing with the vector(s) hold for the vector(s) of weights.
\li Names are a copy of the names for \c m1, with the names for \c m2 appended to the row or column list, as appropriate.
\li This function uses the \ref designated syntax for inputs.
\ingroup data_struct
*/
APOP_VAR_HEAD apop_data *apop_data_stack(apop_data *m1, apop_data * m2, char posn, char inplace){
apop_data * apop_varad_var(m1, NULL)
apop_data * apop_varad_var(m2, NULL)
char apop_varad_var(posn, 'r')
Apop_stopif(!(posn == 'r' || posn == 'c'), return NULL, 0, "Valid positions are 'r' or 'c'"
" you gave me '%c'. Returning NULL.", posn);
char apop_varad_var(inplace, 'n')
inplace = (inplace == 'y' || inplace == 1 || inplace == 'Y') ? 1 : 0;
APOP_VAR_ENDHEAD
if (!m1) return apop_data_copy(m2);
if (!m2) return inplace ? m1 : apop_data_copy(m1);
apop_data *out = NULL;
if (inplace)
out = m1;
else {
apop_data *m = m1->more; //not following the more pointer.
m1->more =NULL;
out = apop_data_copy(m1);
Apop_stopif(out->error, return out, 0, "initial copy failed; leaving.");
m1->more = m;
}
Get_vmsizes(m1); //original sizes of vsize, msize1, msize2.
if (m2->names && !out->names) out->names = apop_name_alloc();
if (posn == 'c'){
if (m2->vector && out->vector){
gsl_matrix_view mview = gsl_matrix_view_vector(m2->vector, m2->vector->size, 1);
out->matrix = apop_matrix_stack(out->matrix, &mview.matrix, posn, .inplace='y');
apop_name_stack(out->names, m2->names, 'c', 'v');
if (m2->names && !m2->names->vector && m2->names->colct) apop_name_add(out->names, "v", 'c');
}
if (m2->vector && !out->vector) {
out->vector= apop_vector_copy(m2->vector);
if (m2->names->vector) apop_name_add(out->names, m2->names->vector, 'v');
}
}
out->matrix = apop_matrix_stack(out->matrix, m2->matrix, posn, .inplace='y');
if (posn == 'r'){
out->vector = apop_vector_stack(out->vector, m2->vector, .inplace='y');
out->weights = apop_vector_stack(out->weights, m2->weights, .inplace='y');
}
if (m2->text){ //we've already copied m1->text, if any, so if m2->text is NULL, we're done.
if (posn=='r'){
Apop_stopif(out->text && m2->textsize[1]!=out->textsize[1],
out->error='d'; return out, 0,
"The first data set has %zu columns of text and the second has %zu columns. "
"I can't stack that.", out->textsize[1], m2->textsize[1]);
int basetextsize = out->textsize[0];
apop_text_alloc(out, basetextsize+m2->textsize[0], out->textsize[1]);
Apop_stopif(out->error, return out, 0, "Allocation error.");
for(int i=0; i< m2->textsize[0]; i++)
for(int j=0; j< m2->textsize[1]; j++)
if (m2->text[i][j] == apop_nul_string)
apop_text_blank(out, i+basetextsize, j);
else apop_text_add(out, i+basetextsize, j, m2->text[i][j]);
} else {
Apop_stopif(out->text && m2->textsize[0]!=out->textsize[0],
out->error='d'; return out, 0,
"The first data set has %zu rows of text and the second has %zu rows. "
"I can't stack that.", out->textsize[0], m2->textsize[0]);
int basetextsize = out->textsize[1];
apop_text_alloc(out, out->textsize[0], basetextsize+m2->textsize[1]);
Apop_stopif(out->error, out->error='a'; return out, 0, "Allocation error.");
for(int i=0; i< m2->textsize[0]; i++)
for(int j=0; j< m2->textsize[1]; j++)
if (m2->text[i][j] == apop_nul_string)
apop_text_blank(out, i, j+basetextsize);
else apop_text_add(out, i, j+basetextsize, m2->text[i][j]);
apop_name_stack(out->names, m2->names, 't');
}
}
if ((posn=='r' && m2->names && m2->names->rowct) || (posn=='c' && m2->names && m2->names->colct)){
int min = posn =='r' ? m1->names->rowct : m1->names->colct;
int max = posn =='r' ? GSL_MAX(vsize, msize1) : msize2;
for (int k = min; k < max; k++) //pad so the name stacking is aligned (if needed)
apop_name_add(out->names, "", posn);
apop_name_stack(out->names, m2->names, posn);
}
return out;
}
/** Split one input \ref apop_data structure into two.
For the opposite operation, see \ref apop_data_stack.
\param in The \ref apop_data structure to split
\param splitpoint The index of what will be the first row/column of the second data set. E.g., if this is -1 and \c 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 \c 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).
\param r_or_c If 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.
\return An array of two \ref apop_data sets. If one is empty then a
\c NULL pointer will be returned in that position. For example, for a data set of 50 rows, <tt>apop_data **out = apop_data_split(data, 100, 'r')</tt> sets <tt>out[0] = apop_data_copy(data)</tt> and <tt>out[1] = NULL</tt>.
\li When splitting at a row, the text is also split.
\li \c more pointer is ignored.
\li The <tt>apop_data->vector</tt> is taken to be the -1st element of the matrix.
\li 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.
\li Data is copied, so you may want to call <tt>apop_data_free(in)</tt> after this.
*/
apop_data ** apop_data_split(apop_data *in, int splitpoint, char r_or_c){
//A long, dull series of contingencies. Bonus: a reasonable use of goto.
apop_data **out = malloc(2*sizeof(apop_data *));
out[0] = out[1] = NULL;
Apop_stopif(!in, return out, 1, "input was NULL; output will be an array of two NULLs.");
gsl_vector v1, v2, w1, w2;
gsl_matrix m1, m2;
int set_v1 = 1, set_v2 = 1,
set_m1 = 1, set_m2 = 1,
set_w1 = 1, set_w2 = 1,
namev0 = 0, namev1 = 0,
namer0 = 0, namer1 = 0,
namec0 = 0, namec1 = 0,
namersplit = -1, namecsplit = -1;
if (r_or_c == 'r' || r_or_c == 'R') {
if (splitpoint <=0)
out[1] = apop_data_copy(in);
else if (in->matrix && splitpoint >= in->matrix->size1)
out[0] = apop_data_copy(in);
else {
namev0 =
namev1 =
namec0 =
namec1 = 1;
if (in->vector){
v1 = gsl_vector_subvector(in->vector, 0, splitpoint).vector;
v2 = gsl_vector_subvector(in->vector, splitpoint, in->vector->size - splitpoint).vector;
} else
set_v1 =
set_v2 = 0;
if (in->weights){
w1 = gsl_vector_subvector(in->weights, 0, splitpoint).vector;
w2 = gsl_vector_subvector(in->weights, splitpoint,
in->weights->size - splitpoint).vector;
} else
set_w1 =
set_w2 = 0;
if (in->matrix){
m1 = gsl_matrix_submatrix (in->matrix, 0, 0, splitpoint, in->matrix->size2).matrix;
m2 = gsl_matrix_submatrix (in->matrix, splitpoint, 0,
in->matrix->size1 - splitpoint, in->matrix->size2).matrix;
} else
set_m1 =
set_m2 = 0;
namersplit=splitpoint;
goto allocation;
}
} else if (r_or_c == 'c' || r_or_c == 'C') {
if (in->weights){
w1 = gsl_vector_subvector(in->weights, 0, in->weights->size).vector;
w2 = gsl_vector_subvector(in->weights, 0, in->weights->size).vector;
} else
set_w1 =
set_w2 = 0;
namer0 = 1;
namer1 = 1;
if (splitpoint <= -1)
out[1] = apop_data_copy(in);
else if (in->matrix && splitpoint >= in->matrix->size2)
out[0] = apop_data_copy(in);
else if (splitpoint == 0){
if (in->vector){
v1 = gsl_vector_subvector(in->vector, 0, in->vector->size).vector;
namev0 = 1;
} else
set_v1 = 0;
set_v2 = 0;
set_m1 = 0;
if (in->matrix){
m2 = gsl_matrix_submatrix (in->matrix, 0, 0,
in->matrix->size1, in->matrix->size2).matrix;
namec1 = 1;
} else
set_m2 = 0;
goto allocation;
} else if (splitpoint > 0 && in->matrix && splitpoint < in->matrix->size2){
if (in->vector){
v1 = gsl_vector_subvector(in->vector, 0, in->vector->size).vector;
namev0 = 1;
} else
set_v1 = 0;
set_v2 = 0;
if (in->matrix){
m1 = gsl_matrix_submatrix (in->matrix, 0, 0, in->matrix->size1, splitpoint).matrix;
m2 = gsl_matrix_submatrix (in->matrix, 0, splitpoint,
in->matrix->size1, in->matrix->size2-splitpoint).matrix;
namecsplit = splitpoint;
} else
set_m1 =
set_m2 = 0;
goto allocation;
} else { //splitpoint >= in->matrix->size2
if (in->vector){
v1 = gsl_vector_subvector(in->vector, 0, in->vector->size).vector;
namev0 = 1;
} else
set_v1 = 0;
set_v2 = 0;
if (in->matrix){
m1 = gsl_matrix_submatrix (in->matrix, 0, 0,
in->matrix->size1, in->matrix->size2).matrix;
namec0 = 1;
}
else set_m1 = 0;
set_m2 = 0;
goto allocation;
}
} else Apop_notify(0, "Please set r_or_c == 'r' or == 'c'. Returning two NULLs.");
return out;
allocation:
out[0] = apop_data_alloc();
out[1] = apop_data_alloc();
if (set_v1) out[0]->vector = apop_vector_copy(&v1);
if (set_v2) out[1]->vector = apop_vector_copy(&v2);
if (set_m1) out[0]->matrix = apop_matrix_copy(&m1);
if (set_m2) out[1]->matrix = apop_matrix_copy(&m2);
if (set_w1) out[0]->weights = apop_vector_copy(&w1);
if (set_w2) out[1]->weights = apop_vector_copy(&w2);
if (namev0 && out[0]) apop_name_stack(out[0]->names, in->names, 'v');
if (namev1 && out[1]) apop_name_stack(out[1]->names, in->names, 'v');
if (namersplit >=0)
for (int k=0; k< in->names->rowct; k++){
int which = (k >= namersplit);
assert(out[which]);
apop_name_add(out[which]->names, in->names->row[k], 'r');
}
else {
if (namer0 && out[0]) apop_name_stack(out[0]->names, in->names, 'r');
if (namer1 && out[1]) apop_name_stack(out[1]->names, in->names, 'r');
}
if (namecsplit >=0)
for (int k=0; k< in->names->colct; k++){
int which = (k >= namecsplit);
assert(out[which]);
apop_name_add(out[which]->names, in->names->col[k], 'c');
}
else {
if (namec0 && out[0]) apop_name_stack(out[0]->names, in->names, 'c');
if (namec1 && out[1]) apop_name_stack(out[1]->names, in->names, 'c');
}
//finally, the text [split by rows only]
if (r_or_c=='r' && in->textsize[0] && in->textsize[1]){
apop_name_stack(out[1]->names, in->names, 't');
apop_text_alloc(out[0], splitpoint, in->textsize[1]);
Apop_stopif(out[0]->error, return out, 0, "Allocation error.");
if (in->textsize[0] > splitpoint){
apop_name_stack(out[0]->names, in->names, 't');
apop_text_alloc(out[1], in->textsize[0]-splitpoint, in->textsize[1]);
Apop_stopif(out[1]->error, return out, 0, "Allocation error.");
}
for (int i=0; i< in->textsize[0]; i++)
for (int j=0; j< in->textsize[1]; j++){
int whichtext = (i >= splitpoint);
int row = whichtext ? i - splitpoint : i;
Asprintf(&(out[whichtext]->text[row][j]), "%s", in->text[i][j]);
}
}
return out;
}
/** Remove the columns set to one in the \c drop vector.
\param n the \ref apop_name structure to be pared down
\param drop a vector with n->colct elements, mostly zero, with a one marking those columns to be removed.
\ingroup names
*/
static void apop_name_rm_columns(apop_name *n, int *drop){
apop_name *newname = apop_name_alloc();
size_t initial_colct = n->colct;
for (size_t i=0; i< initial_colct; i++){
if (drop[i]==0) apop_name_add(newname, n->col[i],'c');
else n->colct--;
free(n->col[i]);
}
free(n->col);
n->col = newname->col;
//we need to free the newname struct, but leave the column intact.
newname->col = NULL;
newname->colct = 0;
apop_name_free(newname);
}
/** Remove the columns set to one in the \c 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.
\param d the \ref apop_data structure to be pared down.
\param drop an array of ints. If use[7]==1, then column seven will be cut from the
output. A reminder: <tt>calloc(in->size2 , sizeof(int))</tt> will fill your array with zeros on allocation, and
<tt>memset(use, 1, in->size2 * sizeof(int))</tt> will
quickly fill an array of ints with nonzero values.
\ingroup data_struct
*/
void apop_data_rm_columns(apop_data *d, int *drop){
gsl_matrix *freeme = d->matrix;
d->matrix = apop_matrix_rm_columns(d->matrix, drop);
gsl_matrix_free(freeme);
apop_name_rm_columns(d->names, drop);
}
/** \def apop_data_prune_columns(in, ...)
Keep only the columns of a data set that you name.
\param in The data set to prune.
\param ... A list of names to retain (i.e. the columns that shouldn't be pruned
out). For example, if you have run \ref apop_data_summarize, you have columns for several
statistics, but may care about only one or two; see the example.
For example:
\include test_pruning.c
\li I use a case-insensitive search to find your column.
\li If your name multiple columns, I'll only give you the first.
\li If I can't find a column matching one of your strings, I throw an error to the screen and continue.
\li This is a macro calling \ref apop_data_prune_columns_base. It packages your list of
columns into a list of strings, adds a \c NULL string at the end, and calls that function.
\hideinitializer */
/** Keep only the columns of a data set that you name.
This is the function called internally by the \ref apop_data_prune_columns macro. In
most cases, you'll want to use that macro. An example of the two uses demonstrating the
difference:
\code
apop_data_prune_columns(d, "mean", "median");
char *list[] = {"mean", "median", NULL};
apop_data_prune_columns_base(d, list);
\endcode
\param d The data set to prune.
\param colnames A null-terminated list of names to retain (i.e. the columns that shouldn't be pruned
out).
\return A pointer to the input data set, now pruned.
*/
apop_data* apop_data_prune_columns_base(apop_data *d, char **colnames){
/* In types.h, you'll find an alias that takes the input, wraps it in the cruft that is
C's compound literal syntax, and appends a final "" to the list of strings. Here, I
find each element of the list, using that "" as a stopper, and then call apop_data_rm_columns.*/
Apop_stopif(!d, return NULL, 1, "You're asking me to prune a NULL data set; returning.");
Apop_stopif(!d->matrix, return d, 1, "You're asking me to prune a data set with NULL matrix; returning.");
int rm_list[d->matrix->size1];
int keep_count = 0;
char **name_step = colnames;
//to throw errors for typos (and slight efficiency gains), I need an array of whether
//each input colname has been used.
while (*name_step++)
keep_count++;
int used_field[keep_count];
memset(used_field, 0, keep_count*sizeof(int));
for (int i=0; i< d->names->colct; i++){
int keep = 0;
for (int j=0; j<keep_count; j++)
if (!used_field[j] && !strcasecmp(d->names->col[i], colnames[j])){
keep ++;
used_field[j]++;
break;
}
rm_list[i] = !keep;
}
apop_data_rm_columns(d, rm_list);
for (int j=0; j<keep_count; j++)
Apop_stopif(!used_field[j], , 1, "You asked me to keep column \"%s\" but I couldn't find a match for it. Typo?", colnames[j]);
return d;
}
/** \defgroup data_set_get Set/get/point to the data element at the given point
\{
First, some examples:
\code
apop_data *d = apop_data_alloc(10, 10, 10);
apop_name_add(d->names, "Zeroth row", 'r');
apop_name_add(d->names, "Zeroth col", 'c');
apop_data_set(d, 8, 0, 27);
assert(apop_data_get(d, 8, .colname="Zeroth") == 27);
double *x = apop_data_ptr(d, .col=7, .rowname="Zeroth");
*x = 270;
assert(apop_data_get(d, 0, 7) == 270);
//apop_data set holding a scalar:
apop_data *s = apop_data_alloc(1);
apop_data_set(s, .val=12);
assert(apop_data_get(s) == 12);
//apop_data set holding a vector:
apop_data *v = apop_data_alloc(12);
for (int i=0; i< 12; i++) apop_data_set(s, i, .val=i*10);
assert(apop_data_get(s,3) == 30);
\endcode
A call like <tt> apop_data_set(in, row, col, data)</tt> is much like the GSL's
<tt>gsl_matrix_set(in->matrix, row, col, data)</tt>,
but with some differences:
\li The \ref apop_data set has names, so we can get/set elements using those names.
\li The versions that take a column/row name use \ref apop_name_find
for the search; see notes there on the name matching rules.
\li The \ref apop_data set has both matrix and vector elements.
\li For those that take a column number, column -1 is the vector element.
\li For those that take a column name, I will search the vector last---if I don't find the name among the matrix columns, but the name matches the vector name, I return column -1.
\li If you give me both a .row and a .rowname, I go with the name; similarly for .col and
.colname.
\li You can give me the name of a page, e.g.
\code
double AIC = apop_data_get(data, .rowname="AIC", .col=-1, .page="<Info>");
\endcode
\li The column (like all defaults) is zero unless stated otherwise, so <tt>apop_data_get(dataset, 1)</tt> gets item (1, 0) from the matrix element of \c dataset. As a do-what-I-mean exception, if there is no matrix element but there is a vector, then this form will get vector element 1. Relying on this DWIM exception is useful iff you can guarantee that a data set will have only a vector or a matrix but not both. Otherwise, be explicit: <tt>apop_data_get(dataset, 1, -1)</tt>.
The \c _ptr functions return a pointer to the given cell. Those functions follow the lead of \c gsl_vector_ptr and \c gsl_matrix_ptr, and like those functions, return a pointer to the appropriate \c double.
\li These functions use the \ref designated syntax for inputs.
*/
/** Get a pointer to an element of an \ref apop_data set.
\li If a \c NULL vector or matrix (as the case may be), stop (unless <tt>apop_opts.stop_on_warning='n'</tt>, then return \c NULL).
\li If the row/column you requested is outside the bounds of the matrix (or the name isn't found), always return \c NULL.
\li See \ref data_set_get "the set/get page" for details.
\param data The data set. Must not be \c NULL.
\param row The row number of the desired element. If <tt>rowname==NULL</tt>, default is zero.
\param col The column number of the desired element. -1 indicates the vector. If <tt>colname==NULL</tt>, default is zero.
\param rowname The row name of the desired element. If <tt>NULL</tt>, use the row number.
\param colname The column name of the desired element. If <tt>NULL</tt>, use the column number.
\param page The case-insensitive name of the page on which the element is found. If \c NULL, use first page.
\return A pointer to the element.
*/
APOP_VAR_HEAD double * apop_data_ptr(apop_data *data, int row, int col, const char *rowname, const char *colname, const char *page){
apop_data * apop_varad_var(data, NULL);
Apop_stopif(!data, return NULL, 0, "You sent me a NULL data set. Returning NULL pointer.");
int apop_varad_var(row, 0);
int apop_varad_var(col, 0);
const char * apop_varad_var(rowname, NULL);
const char * apop_varad_var(colname, NULL);
const char * apop_varad_var(page, NULL);
if (page){
data = apop_data_get_page(data, page);
Apop_stopif(!data, return NULL, 1, "I couldn't find a page with label '%s'. Returning NULL.", page);
};
if (rowname){
row = apop_name_find(data->names, rowname, 'r');
Apop_stopif(row == -2, return NULL, 1, "Couldn't find '%s' amongst the row names.", rowname);
}
if (colname){
col = apop_name_find(data->names, colname, 'c');
Apop_stopif(col == -2, return NULL, 1, "Couldn't find '%s' amongst the column names.", colname);
}
APOP_VAR_ENDHEAD
if (col == -1 || (col == 0 && !data->matrix && data->vector)){
Apop_stopif(!data->vector, return NULL, 1, "You asked for the vector element (col=-1) but it is NULL. Returning NULL.");
return gsl_vector_ptr(data->vector, row);
} else {
Apop_stopif(!data->matrix, return NULL, 1, "You asked for the matrix element (%i, %i) but the matrix is NULL Returning NULL..", row, col);
return gsl_matrix_ptr(data->matrix, row,col);
}
return NULL;//the main function is blank.
}
/** Returns the data element at the given point.
In case of error (probably that you asked for a data point out of bounds), returns \c GSL_NAN.
See \ref data_set_get "the set/get page" for details.
\param data The data set. Must not be \c NULL.
\param row The row number of the desired element. If <tt>rowname==NULL</tt>, default is zero.
\param col The column number of the desired element. -1 indicates the vector. If <tt>colname==NULL</tt>, default is zero.
\param rowname The row name of the desired element. If <tt>NULL</tt>, use the row number.
\param colname The column name of the desired element. If <tt>NULL</tt>, use the column number.
\param page The case-insensitive name of the page on which the element is found. If \c NULL, use first page.
\return The value at the given location. */
APOP_VAR_HEAD double apop_data_get(const apop_data *data, size_t row, int col, const char *rowname, const char *colname, const char *page){
const apop_data * apop_varad_var(data, NULL);
Apop_stopif(!data, return NAN, 0, "You sent me a NULL data set. Returning NaN.");
size_t apop_varad_var(row, 0);
int apop_varad_var(col, 0);
const char * apop_varad_var(rowname, NULL);
const char * apop_varad_var(colname, NULL);
const char * apop_varad_var(page, NULL);
if (page){
data = apop_data_get_page(data, page);
Apop_stopif(!data, return NAN, 1, "I couldn't find a page with label '%s'. Returning NaN.", page);
};
if (rowname){
row = apop_name_find(data->names, rowname, 'r');
Apop_stopif(row == -2, return NAN, 1, "Couldn't find '%s' amongst the row names. Returning NaN.", rowname);
}
if (colname){
col = apop_name_find(data->names, colname, 'c');
Apop_stopif(col == -2, return NAN, 1, "Couldn't find '%s' amongst the column names. Returning NaN.", colname);
}
APOP_VAR_ENDHEAD
if (col==-1 || (col == 0 && !data->matrix && data->vector)){
Apop_stopif(!data->vector, return NAN, 1, "You asked for the vector element (col=-1) but it is NULL.");
return gsl_vector_get(data->vector, row);
} else {
Apop_stopif(!data->matrix, return NAN, 1, "You asked for the matrix element (%zu, %i) but the matrix is NULL.", row, col);
return gsl_matrix_get(data->matrix, row, col);
}
}
/* The only hint the GSL gives that something failed is that the error-handler is called.
The error handling function won't let you set an output to the function. So all we
can do is use a global variable.
*/
static threadlocal int error_for_set; //see apop_internal.h
void apop_gsl_error_for_set(const char *reason, const char *file, int line, int gsl_errno){
Apop_notify(1, "%s: %s", file, reason);
Apop_maybe_abort(1);
error_for_set = -1;
}
/** Set a data element.
This function uses the \ref designated syntax, so with names you don't have to worry
about element ordering. But the ordering of elements may still be
noteworthy. For compatibility with older editions of Apophenia, the order is (row, col,
value, rowname, colname), so the following would all set row 3, column 8, of \c d to 5:
\code
apop_data_set(d, 3, 8, 5);
apop_data_set(d, .row = 3, .col=8, .val=5);
apop_data_set(d, .row = 3, .colname="Column 8", .val=5);
//but:
apop_data_set(d, .row = 3, .colname="Column 8", 5); //invalid---the value doesn't follow the colname.
\endcode
\return 0=OK, -1=error (couldn't find row/column name, or you asked for a location outside the vector/matrix bounds).
\li The error codes for out-of-bounds errors are thread-safe iff you are have a
C11-compliant compiler (thanks to the \c _Thread_local keyword) or a version of GCC with the \c __thread
extension enabled.
See \ref data_set_get "the set/get page" for details.
\param data The data set. Must not be \c NULL.
\param row The row number of the desired element. If <tt>rowname==NULL</tt>, default is zero.
\param col The column number of the desired element. -1 indicates the vector. If <tt>colname==NULL</tt>, default is zero.
\param rowname The row name of the desired element. If <tt>NULL</tt>, use the row number.
\param colname The column name of the desired element. If <tt>NULL</tt>, use the column number.
\param page The case-insensitive name of the page on which the element is found. If \c NULL, use first page.
\param val The value to give the point.
\return The value at the given location. */
APOP_VAR_HEAD int apop_data_set(apop_data *data, size_t row, int col, const double val, const char *colname, const char *rowname, const char *page){
apop_data * apop_varad_var(data, NULL);
Apop_stopif(!data, return -1, 0, "You sent me a NULL data set.");
size_t apop_varad_var(row, 0);
int apop_varad_var(col, 0);
const double apop_varad_var(val, 0);
const char * apop_varad_var(rowname, NULL);
const char * apop_varad_var(colname, NULL);
const char * apop_varad_var(page, NULL);
if (page){
data = apop_data_get_page((apop_data*)data, page);
Apop_stopif(!data, return -1, 1, "I couldn't find a page with label '%s'. Making no changes.", page);
}
if (rowname){
row = apop_name_find(data->names, rowname, 'r');
Apop_stopif(row == -2, return -1, 1, "Couldn't find '%s' amongst the column names. Making no changes.", rowname);
}
if (colname){
col = apop_name_find(data->names, colname, 'c');
Apop_stopif(col == -2, return -1, 1, "Couldn't find '%s' amongst the column names. Making no changes.", colname);
}
APOP_VAR_ENDHEAD
Set_gsl_handler
if (col==-1 || (col == 0 && !data->matrix && data->vector)){
Apop_stopif(!data->vector, return -1, 1, "You're trying to set a vector element (row=-1) but the vector is NULL.");
gsl_vector_set(data->vector, row, val);
} else {
Apop_stopif(!data->matrix, return -1, 1, "You're trying to set the matrix element (%zu, %i) but the matrix is NULL.", row, col);
gsl_matrix_set(data->matrix, row, col, val);
}
Unset_gsl_handler
return error_for_set;
}
/** \} //End data_set_get group */
/** Now that you've used \ref Apop_r to pull a row from an \ref 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:
\li a vector element has to be present if one existed in the original,
\li same for the weights vector,
\li the matrix in the destination has to have as many columns as in the original, and
\li the text has to have a row long enough to hold the original
\li If the row to be written to already has a rowname, it is overwritten.
If <tt>d->names->rowct == row_number</tt> (all rows up to \c 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.
\li Column names (of all types) aren't touched. Maybe use \c apop_data_copy or \c apop_name_copy if you need to copy these names.
If any of the source elements are \c NULL, I won't bother to check that element in the
destination.
\return 0=OK, -1=error (probably a source/destination size mismatch).
\li The error codes for out-of-bounds errors are thread-safe iff you are have a
C11-compliant compiler (thanks to the \c _Thread_local keyword) or a version of GCC with the \c __thread extension enabled.
*/
int apop_data_set_row(apop_data * d, apop_data *row, int row_number){
Set_gsl_handler
if (row->vector){
Apop_assert_negone(d->vector, "You asked me to copy an apop_data row with a vector element to "
"an apop_data set with no vector.");
gsl_vector_set(d->vector, row_number, row->vector->data[0]);
}
if (row->matrix && row->matrix->size2 > 0){
Apop_assert_negone(d->matrix, "You asked me to copy an apop_data row with a matrix row to "
"an apop_data set with no matrix.");
gsl_vector_memcpy(Apop_rv(d, row_number), Apop_rv(row, 0));
}
if (row->textsize[1]){
Apop_assert_negone(d->textsize[1], "You asked me to copy an apop_data row with text to "
"an apop_data set with no text element.");
for (int i=0; i < row->textsize[1]; i++){
free(d->text[row_number][i]);
d->text[row_number][i]=
row->text[0][i] == apop_nul_string
? apop_nul_string
: strdup(row->text[0][i]);
}
}
if (row->weights){
Apop_assert_negone(d->weights, "You asked me to copy an apop_data row with a weight to "
"an apop_data set with no weights vector.");
gsl_vector_set(d->weights, row_number, row->weights->data[0]);
}
if (row->names && row->names->rowct && d->names){
if (row_number < d->names->rowct){
free(d->names->row[row_number]);
d->names->row[row_number]=strdup(row->names->row[0]);
} else if (row_number == d->names->rowct)
apop_name_add(d->names, row->names->row[0], 'r');
}
Unset_gsl_handler
return error_for_set;
}
/** 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.