Gitlab will migrate to a new storage backend starting 0300 UTC on 2020-04-04. We do not anticipate a maintenance window for this migration. Performance may be impacted over the weekend. Thanks for your patience.

LALInferenceGenerateROQ.c 77.6 KB
Newer Older
1 2 3
/*
 *  LALInferenceCreateROQ.c: Reduced order quadrature basis and interpolant generation
 *
4
 *  Copyright (C) 2014, 2016 Matthew Pitkin, Rory Smith
5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23
 *
 *  This program is free software; you can redistribute it and/or modify
 *  it under the terms of the GNU General Public License as published by
 *  the Free Software Foundation; either version 2 of the License, or
 *  (at your option) any later version.
 *
 *  This program is distributed in the hope that it will be useful,
 *  but WITHOUT ANY WARRANTY; without even the implied warranty of
 *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 *  GNU General Public License for more details.
 *
 *  You should have received a copy of the GNU General Public License
 *  along with with program; see the file COPYING. If not, write to the
 *  Free Software Foundation, Inc., 59 Temple Place, Suite 330, Boston,
 *  MA  02111-1307  USA
 */

#include <lal/LALInferenceGenerateROQ.h>

24 25 26 27 28
#ifndef _OPENMP
#define omp ignore
#endif


29 30 31 32 33 34 35
/* internal function definitions */

/* define function to project model vector onto the training set of models */
void project_onto_basis(gsl_vector *weight,
                        gsl_matrix *RB,
                        gsl_matrix *TS,
                        gsl_matrix *projections,
36
                        INT4 idx);
37 38 39 40 41

void complex_project_onto_basis(gsl_vector *weight,
                                gsl_matrix_complex *RB,
                                gsl_matrix_complex *TS,
                                gsl_matrix_complex *projections,
42
                                INT4 idx);
43 44

/* the dot product of two real vectors scaled by a weighting factor */
45
REAL8 weighted_dot_product(const gsl_vector *weight, gsl_vector *a, gsl_vector *b);
46 47

/* the dot product of two complex vectors scaled by a weighting factor */
48
gsl_complex complex_weighted_dot_product(const gsl_vector *weight, gsl_vector_complex *a, gsl_vector_complex *b);
49

50 51
void normalise(const gsl_vector *weight, gsl_vector *a);
void complex_normalise(const gsl_vector *weight, gsl_vector_complex *a);
52 53 54 55

void normalise_training_set(gsl_vector *weight, gsl_matrix *TS);
void complex_normalise_training_set(gsl_vector *weight, gsl_matrix_complex *TS);

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
double normalisation(const gsl_vector *weight, gsl_vector *a);
double complex_normalisation(const gsl_vector *weight, gsl_vector_complex *a);

void modified_gs_complex(gsl_vector_complex *ru,
                 gsl_vector_complex *ortho_basis,
                 const gsl_matrix_complex *RB_space,
                 const gsl_vector *wQuad,
                 const int dim_RB);

void iterated_modified_gm_complex(gsl_vector_complex *ru,
                                  gsl_vector_complex *ortho_basis,
                                  const gsl_matrix_complex *RB_space,
                                  const gsl_vector *wQuad,
                                  const int dim_RB);

void modified_gs(gsl_vector *ru,
                 gsl_vector *ortho_basis,
                 const gsl_matrix *RB_space,
                 const gsl_vector *wQuad,
                 const int dim_RB);

void iterated_modified_gm(gsl_vector *ru,
                          gsl_vector *ortho_basis,
                          const gsl_matrix *RB_space,
                          const gsl_vector *wQuad,
                          const int dim_RB);

83
/* get the B_matrix */
84 85
REAL8Array *B_matrix(gsl_matrix *V, gsl_matrix *RB);
COMPLEX16Array *complex_B_matrix(gsl_matrix_complex *V, gsl_matrix_complex *RB);
86 87 88 89

/* find the index of the absolute maximum value for a complex vector */
int complex_vector_maxabs_index( gsl_vector_complex *c );

90

91 92 93 94 95 96 97 98 99 100 101 102 103 104
/** \brief Function to project the training set onto a given basis vector
 *
 * This is an internal function to be used by \c LALInferenceCreateREAL8OrthonormalBasis
 *
 * @param[in] weight The normalisation weight(s) for the training set waveforms (e.g. time or frequency step(s) size)
 * @param[in] RB The reduced basis set
 * @param[in] TS The training set of waveforms
 * @param[in] projections The set of projections (this is updated in this function)
 * @param[in] idx The index of the reduced basis vector that the training set will project onto
 */
void project_onto_basis(gsl_vector *weight,
                        gsl_matrix *RB,
                        gsl_matrix *TS,
                        gsl_matrix *projections,
105
                        INT4 idx){
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
  size_t row = 0;
  gsl_vector_view basis;

  XLAL_CALLGSL( basis = gsl_matrix_row(RB, idx) );

  for ( row=0; row < TS->size1; row++ ){
    double prod;
    gsl_vector_view proj, model;
    gsl_vector *basisscale;

    XLAL_CALLGSL( proj = gsl_matrix_row(projections, row) );
    XLAL_CALLGSL( basisscale = gsl_vector_calloc(TS->size2) );

    XLAL_CALLGSL( model = gsl_matrix_row(TS, row) ); /* get model from training set */

    prod = weighted_dot_product(weight, &basis.vector, &model.vector);

    XLAL_CALLGSL( gsl_vector_memcpy(basisscale, &basis.vector) );
    XLAL_CALLGSL( gsl_vector_scale(basisscale, prod) );
    XLAL_CALLGSL( gsl_vector_add(&proj.vector, basisscale) );
    XLAL_CALLGSL( gsl_vector_free(basisscale) );
  }
}


/** \brief Function to project the complex training set onto a given basis vector
 *
 * This is an internal function to be used by \c LALInferenceCreateCOMPLEX16OrthonormalBasis
 *
 * @param[in] weight The normalisation weight(s) for the training set waveforms (e.g. time or frequency step(s) size)
 * @param[in] RB The reduced basis set
 * @param[in] TS The training set of waveforms
 * @param[in] projections The set of projections (this is updated in this function)
 * @param[in] idx The index of the reduced basis vector that the training set will project onto
 */
void complex_project_onto_basis(gsl_vector *weight,
                                gsl_matrix_complex *RB,
                                gsl_matrix_complex *TS,
                                gsl_matrix_complex *projections,
145
                                INT4 idx){
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
  size_t row = 0;
  gsl_vector_complex_view basis;

  XLAL_CALLGSL( basis = gsl_matrix_complex_row(RB, idx) );

  for ( row=0; row < TS->size1; row++ ){
    gsl_complex cprod;
    gsl_vector_complex_view proj, model;
    gsl_vector_complex *basisscale;

    XLAL_CALLGSL( proj = gsl_matrix_complex_row(projections, row) );
    XLAL_CALLGSL( basisscale = gsl_vector_complex_calloc(TS->size2) );

    XLAL_CALLGSL( model = gsl_matrix_complex_row(TS, row) ); /* get model from training set */

    cprod = complex_weighted_dot_product(weight, &basis.vector, &model.vector);

    XLAL_CALLGSL( gsl_vector_complex_memcpy(basisscale, &basis.vector) );
    XLAL_CALLGSL( gsl_vector_complex_scale(basisscale, cprod) );
    XLAL_CALLGSL( gsl_vector_complex_add(&proj.vector, basisscale) );
    XLAL_CALLGSL( gsl_vector_complex_free(basisscale) );
  }
}


/** \brief The dot product of two real vectors scaled by a given weight factor
 *
 * @param[in] weight A (set of) scaling factor(s) for the dot product
 * @param[in] a The first vector
 * @param[in] b The second vector
 *
 * @return The real dot product of the two vectors
 */
179
REAL8 weighted_dot_product(const gsl_vector *weight, gsl_vector *a, gsl_vector *b){
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
  REAL8 dp;
  gsl_vector *weighted;

  XLAL_CHECK_REAL8( a->size == b->size, XLAL_EFUNC, "Size of input vectors are not the same.");

  XLAL_CALLGSL( weighted = gsl_vector_calloc(a->size) );
  XLAL_CALLGSL( gsl_vector_memcpy(weighted, a) );

  /* multiply vector by weight */
  if ( weight->size == 1 ){ /* just a single weight to scale with */
    XLAL_CALLGSL( gsl_vector_scale(weighted, gsl_vector_get(weight, 0)) );
  }
  else if ( weight->size == a->size ){
    XLAL_CALLGSL( gsl_vector_mul(weighted, weight) );
  }
  else{
    XLAL_ERROR_REAL8( XLAL_EFUNC, "Vector of weights must either contain a single value, or be the same length as the other input vectors." );
  }

  /* get dot product */
  XLAL_CALLGSL( gsl_blas_ddot(weighted, b, &dp) );

  XLAL_CALLGSL( gsl_vector_free(weighted) );

  return dp;
}


/** \brief The dot product of two complex vectors scaled by a given weight factor
 *
 * The dot product is produced using the complex conjugate of the first vector.
 *
 * @param[in] weight A real scaling factor for the dot product
 * @param[in] a The first complex vector
 * @param[in] b The second complex vector
 *
 * @return The absolute value of the complex dot product of the two vectors
 */
218
gsl_complex complex_weighted_dot_product(const gsl_vector *weight, gsl_vector_complex *a, gsl_vector_complex *b){
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
  gsl_complex dp;
  gsl_vector_complex *weighted;

  if ( a->size != b->size ){ XLAL_PRINT_ERROR( "Size of input vectors are not the same." ); }

  XLAL_CALLGSL( weighted = gsl_vector_complex_calloc(a->size) );
  XLAL_CALLGSL( gsl_vector_complex_memcpy(weighted, a) );

  /* multiply vector by weight */
  if ( weight->size == 1 ){ /* just a single weight to scale with */
    XLAL_CALLGSL( gsl_blas_zdscal(gsl_vector_get(weight, 0), weighted) );
  }
  else if ( weight->size == a->size ){
    gsl_vector_view rview, iview;

    XLAL_CALLGSL( rview = gsl_vector_complex_real(weighted) );
    XLAL_CALLGSL( iview = gsl_vector_complex_imag(weighted) );

    XLAL_CALLGSL( gsl_vector_mul(&rview.vector, weight) );
    XLAL_CALLGSL( gsl_vector_mul(&iview.vector, weight) );
  }
  else{
    XLAL_PRINT_ERROR( "Vector of weights must either contain a single value, or be the same length as the other input vectors." );
  }

  /* get dot product (taking the complex conjugate of the first vector) */
  XLAL_CALLGSL( gsl_blas_zdotc(weighted, b, &dp) );
  XLAL_CALLGSL( gsl_vector_complex_free(weighted) );

  return dp;
}


/** \brief Normalise a real vector with a given weighting
 *
 * @param[in] weight The weighting(s) in the normalisation (e.g. time of frequency step(s) between points)
 * @param[in] a The vector to be normalise (this will be change by the function to return the
 * normalised vector.
 */
258
void normalise(const gsl_vector *weight, gsl_vector *a){
259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277
  double norm;

  if ( weight->size == 1 ){
    XLAL_CALLGSL( norm = gsl_blas_dnrm2(a) ); /* use GSL normalisation calculation function */
    XLAL_CALLGSL( gsl_vector_scale(a, 1./(norm*sqrt(gsl_vector_get(weight, 0)))) );
  }
  else if ( weight->size == a->size ){
    norm = 1./sqrt(weighted_dot_product(weight, a, a));
    XLAL_CALLGSL( gsl_vector_scale(a, norm) );
  }
  else{
    XLAL_ERROR_VOID( XLAL_EFUNC, "Vector of weights must either contain a single value, or be the same length as the other input vectors." );
  }
}


/** \brief Normalise a complex vector with a given (real) weighting
 *
 * @param[in] weight The weighting(s) in the normalisation (e.g. time of frequency step(s) between points)
278 279
 * @param[in] a The vector to be normalised (this will be change by the function to return the
 * normalised vector).
280
 */
281
void complex_normalise(const gsl_vector *weight, gsl_vector_complex *a){
282 283 284 285 286 287 288 289 290 291 292 293 294 295 296
  double norm;

  if ( weight->size == 1 ){
    XLAL_CALLGSL( norm = gsl_blas_dznrm2(a) ); /* use GSL normalisation calculation function */
    XLAL_CALLGSL( gsl_blas_zdscal(1./(norm*sqrt(gsl_vector_get(weight, 0))), a) );
  }
  else if ( weight->size == a->size ){
    norm = 1./sqrt(gsl_complex_abs(complex_weighted_dot_product(weight, a, a)));
    XLAL_CALLGSL( gsl_blas_zdscal(norm, a) );
  }
  else{
    XLAL_ERROR_VOID( XLAL_EFUNC, "Vector of weights must either contain a single value, or be the same length as the other input vectors." );
  }
}

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
/** \brief Get normalisation constant for a real vector
 *
 * @param[in] weight The weighting(s) in the normalisation (e.g. time of frequency step(s) between points)
 * @param[in] a The vector to calculate the normalisation constant for
 *
 * @return The normalisation
 */
double normalisation(const gsl_vector *weight, gsl_vector *a){
  double norm;

  if ( weight->size == 1 ){
    XLAL_CALLGSL( norm = gsl_blas_dnrm2(a) ); /* use GSL normalisation calculation function */
    norm *= sqrt(gsl_vector_get(weight, 0));
  }
  else if ( weight->size == a->size ){
    norm = sqrt(fabs(weighted_dot_product(weight, a, a)));
  }
  else{
    XLAL_ERROR_REAL8( XLAL_EFUNC, "Vector of weights must either contain a single value, or be the same length as the other input vectors." );
  }

  return norm;
}


/** \brief Get normalisation constant for a complex vector
 *
 * @param[in] weight The weighting(s) in the normalisation (e.g. time of frequency step(s) between points)
 * @param[in] a The vector to calculate the normalisation constant for
 *
 * @return The normalisation
 */
double complex_normalisation(const gsl_vector *weight, gsl_vector_complex *a){
  double norm;

  if ( weight->size == 1 ){
    XLAL_CALLGSL( norm = gsl_blas_dznrm2(a) ); /* use GSL normalisation calculation function */
    norm *= sqrt(gsl_vector_get(weight, 0));
  }
  else if ( weight->size == a->size ){
    norm = sqrt(gsl_complex_abs(complex_weighted_dot_product(weight, a, a)));
  }
  else{
    XLAL_ERROR_REAL8( XLAL_EFUNC, "Vector of weights must either contain a single value, or be the same length as the other input vectors." );
  }

  return norm;
}


348 349 350
/** \brief Normalise the set of training waveforms
 *
 * This function will normalise a set of training waveforms. This will be used within the
351
 * \c LALInferenceCreateREAL8OrthonormalBasis function.
352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368
 *
 * @param[in] weight The e.g. time/frequency step in the waveforms used to normalise the waveforms
 * @param[in] TS The training set to be normalised.
 */
void normalise_training_set(gsl_vector *weight, gsl_matrix *TS){
  gsl_vector_view rowview;
  size_t i = 0;
  for ( i=0; i<TS->size1; i++ ){
    XLAL_CALLGSL( rowview = gsl_matrix_row(TS, i) );
    normalise(weight, &rowview.vector);
  }
}


/** \brief Normalise the set of complex training waveforms
 *
 * This function will normalise a set of complex training waveforms. This will be used within the
369
 * \c LALInferenceCreateCOMPLEX16OrthonormalBasis function.
370 371 372 373 374 375 376
 *
 * @param[in] weight The e.g. time/frequency step in the waveforms used to normalise the waveforms
 * @param[in] TS The training set to be normalised.
 */
void complex_normalise_training_set(gsl_vector *weight, gsl_matrix_complex *TS){
  gsl_vector_complex_view rowview;
  size_t i = 0;
377

378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394
  for ( i=0; i<TS->size1; i++ ){
    XLAL_CALLGSL( rowview = gsl_matrix_complex_row(TS, i) );
    complex_normalise(weight, &rowview.vector);
  }
}


/** \brief Get the interpolant of a reduced basis set
 *
 * This is used internally by \c LALInferenceCreateREALROQInterpolant when
 * iteratively calculating the interpolation matrix.
 *
 * @param[in] V The matrix containing the basis vector points at the current interpolation nodes
 * @param[in] RB The set of basis vectors
 *
 * @return The interpolant matrix
 */
395
REAL8Array *B_matrix(gsl_matrix *V, gsl_matrix *RB){
396 397
  /* get inverse of V */
  size_t n = V->size1;
398
  gsl_matrix *invV, *LU;
399 400
  gsl_permutation *p;
  gsl_matrix_view subRB;
401 402
  REAL8Array *B = NULL;
  UINT4Vector *dims = NULL;
403 404 405 406 407 408 409 410 411 412 413 414 415 416 417

  XLAL_CALLGSL( invV = gsl_matrix_alloc(n, n) );
  int signum;

  /* use LU decomposition to get inverse */
  XLAL_CALLGSL( LU = gsl_matrix_alloc(n, n) );
  XLAL_CALLGSL( gsl_matrix_memcpy(LU, V) );

  XLAL_CALLGSL( p = gsl_permutation_alloc(n) );
  XLAL_CALLGSL( gsl_linalg_LU_decomp(LU, p, &signum) );
  XLAL_CALLGSL( gsl_linalg_LU_invert(LU, p, invV) );
  XLAL_CALLGSL( gsl_permutation_free(p) );
  XLAL_CALLGSL( gsl_matrix_free(LU) );

  /* get B matrix */
418 419 420 421 422 423 424
  dims = XLALCreateUINT4Vector( 2 );
  dims->data[0] = n;
  dims->data[1] = RB->size2;
  B = XLALCreateREAL8Array( dims );
  XLALDestroyUINT4Vector( dims );
  gsl_matrix_view Bview;
  Bview = gsl_matrix_view_array(B->data, n, RB->size2);
425
  XLAL_CALLGSL( subRB = gsl_matrix_submatrix(RB, 0, 0, n, RB->size2) );
426
  XLAL_CALLGSL( gsl_blas_dgemm(CblasTrans, CblasNoTrans, 1.0, invV, &subRB.matrix, 0., &Bview.matrix) );
427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443

  XLAL_CALLGSL( gsl_matrix_free(invV) );

  return B;
}


/** \brief Get the interpolant of a complex reduced basis set
 *
 * This is used internally by \c LALInferenceCreateCOMPLEXROQInterpolant when
 * iteratively calculating the interpolation matrix.
 *
 * @param[in] V The matrix containing the basis vector points at the current interpolation nodes
 * @param[in] RB The set of basis vectors
 *
 * @return The interpolant matrix
 */
444
COMPLEX16Array *complex_B_matrix(gsl_matrix_complex *V, gsl_matrix_complex *RB){
445 446
  /* get inverse of V */
  size_t n = V->size1;
447
  gsl_matrix_complex *invV, *LU;
448 449
  gsl_permutation *p;
  gsl_matrix_complex_view subRB;
450 451
  COMPLEX16Array *B = NULL;
  UINT4Vector *dims = NULL;
452 453 454 455 456 457 458 459 460 461 462 463 464 465 466

  XLAL_CALLGSL( invV = gsl_matrix_complex_alloc(n, n) );
  int signum;

  /* use LU decomposition to get inverse */
  XLAL_CALLGSL( LU = gsl_matrix_complex_alloc(n, n) );
  XLAL_CALLGSL( gsl_matrix_complex_memcpy(LU, V) );

  XLAL_CALLGSL( p = gsl_permutation_alloc(n) );
  XLAL_CALLGSL( gsl_linalg_complex_LU_decomp(LU, p, &signum) );
  XLAL_CALLGSL( gsl_linalg_complex_LU_invert(LU, p, invV) );
  XLAL_CALLGSL( gsl_permutation_free(p) );
  XLAL_CALLGSL( gsl_matrix_complex_free(LU) );

  /* get B matrix */
467 468 469 470 471 472 473
  dims = XLALCreateUINT4Vector( 2 );
  dims->data[0] = n;
  dims->data[1] = RB->size2;
  B = XLALCreateCOMPLEX16Array( dims );
  XLALDestroyUINT4Vector( dims );
  gsl_matrix_complex_view Bview;
  Bview = gsl_matrix_complex_view_array((double *)B->data, n, RB->size2);
474
  XLAL_CALLGSL( subRB = gsl_matrix_complex_submatrix(RB, 0, 0, n, RB->size2) );
475
  XLAL_CALLGSL( gsl_blas_zgemm(CblasTrans, CblasNoTrans, GSL_COMPLEX_ONE, invV, &subRB.matrix, GSL_COMPLEX_ZERO, &Bview.matrix) );
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

  XLAL_CALLGSL( gsl_matrix_complex_free(invV) );

  return B;
}

/** \brief Get the index of the maximum absolute value for a complex vector
 *
 * @param[in] c A complex vector
 *
 * @return The index of the maximum absolute value of that vector
 */
int complex_vector_maxabs_index( gsl_vector_complex *c ){
  double maxv = -INFINITY, absval = 0.;
  int idx = 0;
  size_t i = 0;

  for ( i=0; i<c->size; i++ ){
    XLAL_CALLGSL( absval = gsl_complex_abs(gsl_vector_complex_get(c, i)) );

    if ( absval > maxv ){
      maxv = absval;
      idx = (int)i;
    }
  }

  return idx;
}

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
/** \brief Modified Gram-Schmidt algorithm for complex data
 *
 * A modified Gram-Schmidt algorithm taken from the
 * <a href="https://bitbucket.org/sfield83/greedycpp">greedycpp</a> code.
 *
 * @param[out] ru A 1-by-(\c dim_RB + 1) complex vector returning a slice of the R matrix 
 * @param[out] ortho_basis A complex vector returning the orthogonal basis vector
 * @param[in] RB_space A complex matrix containing the existing set of orthonormal basis vectors
 * @param[in] wQuad A vector of weights for the inner products
 * @param[in] dim_RB The number of elements in the current \c RB_space
 */
void modified_gs_complex(gsl_vector_complex *ru,
                         gsl_vector_complex *ortho_basis,
                         const gsl_matrix_complex *RB_space,
                         const gsl_vector *wQuad,
                         const int dim_RB){
  INT4 quad_num = RB_space->size2;
  gsl_complex L2_proj;
  gsl_vector_complex *basis;

  basis = gsl_vector_complex_alloc(quad_num);

  /* if not done, R matrix fills up below diagonal with .5 instead of 0 */
  gsl_vector_complex_set_zero(ru); 

  for(INT4 i = 0; i < dim_RB; i++){
    gsl_matrix_complex_get_row(basis, RB_space, i);

    /* ortho_basis = ortho_basis - L2_proj*basis; */
    L2_proj = complex_weighted_dot_product(wQuad, basis, ortho_basis);
    gsl_vector_complex_set(ru, i, L2_proj);
    gsl_vector_complex_scale(basis, L2_proj); /* basis <- basis*L2_proj */
    gsl_vector_complex_sub(ortho_basis, basis); /* ortho <- ortho - basis */
  }

  double nrm = complex_normalisation(wQuad, ortho_basis);
  gsl_complex nrmc;
  GSL_SET_COMPLEX(&nrmc, nrm, 0.0);
  gsl_vector_complex_set(ru, dim_RB, nrmc);

  complex_normalise(wQuad, ortho_basis);

  gsl_vector_complex_free(basis);
}


/** \brief Modified Gram-Schmidt algorithm for real data
 *
 * A modified Gram-Schmidt algorithm taken from the
 * <a href="https://bitbucket.org/sfield83/greedycpp">greedycpp</a> code.
 *
 * @param[out] ru A 1-by-(\c dim_RB + 1) vector returning a slice of the R matrix 
 * @param[out] ortho_basis A vector returning the orthogonal basis vector
 * @param[in] RB_space A matrix containing the existing set of orthonormal basis vectors
 * @param[in] wQuad A vector of weights for the inner products
 * @param[in] dim_RB The number of elements in the current \c RB_space
 */
void modified_gs(gsl_vector *ru,
                 gsl_vector *ortho_basis,
                 const gsl_matrix *RB_space,
                 const gsl_vector *wQuad,
                 const int dim_RB){
  INT4 quad_num = RB_space->size2;
  REAL8 L2_proj;
  gsl_vector *basis;

  basis = gsl_vector_alloc(quad_num);

  /* if not done, R matrix fills up below diagonal with .5 instead of 0 */
  gsl_vector_set_zero(ru); 

  for(INT4 i = 0; i < dim_RB; i++){
    gsl_matrix_get_row(basis, RB_space, i);

    /* ortho_basis = ortho_basis - L2_proj*basis; */
    L2_proj = weighted_dot_product(wQuad, basis, ortho_basis);
    gsl_vector_set(ru, i, L2_proj);
    gsl_vector_scale(basis, L2_proj); /* basis <- basis*L2_proj */
    gsl_vector_sub(ortho_basis, basis); /* ortho <- ortho - basis */
  }

  double nrm = normalisation(wQuad, ortho_basis);
  gsl_vector_set(ru, dim_RB, nrm);

  normalise(wQuad, ortho_basis);

  gsl_vector_free(basis);
}


/** \brief Iterated modified Gram-Schmidt algorithm for complex data
 *
 * An iterated modified Gram-Schmidt algorithm taken from the
 * <a href="https://bitbucket.org/sfield83/greedycpp">greedycpp</a> code. This uses the method
 * given in \cite Hoffmann1989
 *
 * @param[out] ru A complex 1-by-(\c dim_RB + 1) vector returning a slice of the R matrix 
 * @param[out] ortho_basis A complex vector returning the orthogonal basis vector
 * @param[in] RB_space A complex matrix containing the existing set of orthonormal basis vectors
 * @param[in] wQuad A vector of weights for the inner products
 * @param[in] dim_RB The number of elements in the current \c RB_space
 */
void iterated_modified_gm_complex(gsl_vector_complex *ru,
                                  gsl_vector_complex *ortho_basis,
                                  const gsl_matrix_complex *RB_space,
                                  const gsl_vector *wQuad,
                                  const int dim_RB){
  REAL8 ortho_condition = .5; /* hard coded IMGS stopping condition */

  INT4 quad_num = RB_space->size2;
  INT4 r_size = ru->size;
  REAL8 nrm_prev = complex_normalisation(wQuad, ortho_basis);
  UINT4 flag = 0, iter = 0;
  gsl_vector_complex *e,*r_last;
  REAL8 nrm_current = 0.;
  gsl_complex nrmc_current;

  /* allocate memory */
  e = gsl_vector_complex_alloc(quad_num);
  r_last = gsl_vector_complex_alloc(r_size);

  gsl_vector_complex_memcpy(e,ortho_basis);
  complex_normalise(wQuad, e);

  /* if not done, R matrix fills up below diagonal with .5 instead of 0 */
  gsl_vector_complex_set_zero(ru);  

  while( !flag ){
    gsl_vector_complex_memcpy(ortho_basis, e);
    gsl_vector_complex_set_zero(r_last);

    modified_gs_complex(r_last, ortho_basis, RB_space, wQuad, dim_RB);

    gsl_vector_complex_add(ru, r_last);
    nrmc_current = gsl_vector_complex_get(r_last, dim_RB);
    nrm_current = GSL_REAL(nrmc_current);

    gsl_vector_complex_scale(ortho_basis, nrmc_current);

    if( nrm_current/nrm_prev <= ortho_condition ) {
      nrm_prev = nrm_current;
      iter = iter + 1;
      gsl_vector_complex_memcpy(e, ortho_basis);
    }
    else{ flag = 1; }

    nrm_current = complex_normalisation(wQuad, ortho_basis);
    GSL_SET_COMPLEX(&nrmc_current, nrm_current, 0.0);
    gsl_vector_complex_set(ru, dim_RB, nrmc_current);

    complex_normalise(wQuad, ortho_basis);
  }

  gsl_vector_complex_free(e);
  gsl_vector_complex_free(r_last);
}


/** \brief Iterated modified Gram-Schmidt algorithm for real data
 *
 * An iterated modified Gram-Schmidt algorithm taken from the
 * <a href="https://bitbucket.org/sfield83/greedycpp">greedycpp</a> code. This uses the method
 * given in \cite Hoffmann1989
 *
 * @param[out] ru A 1-by-(\c dim_RB + 1) vector returning a slice of the R matrix 
 * @param[out] ortho_basis A vector returning the orthogonal basis vector
 * @param[in] RB_space A matrix containing the existing set of orthonormal basis vectors
 * @param[in] wQuad A vector of weights for the inner products
 * @param[in] dim_RB The number of elements in the current \c RB_space
 */
void iterated_modified_gm(gsl_vector *ru,
                          gsl_vector *ortho_basis,
                          const gsl_matrix *RB_space,
                          const gsl_vector *wQuad,
                          const int dim_RB){
  REAL8 ortho_condition = .5; /* hard coded IMGS stopping condition */

  INT4 quad_num = RB_space->size2;
  INT4 r_size = ru->size;
  REAL8 nrm_prev = normalisation(wQuad, ortho_basis);
  UINT4 flag = 0, iter = 0;
  gsl_vector *e,*r_last;
  REAL8 nrm_current = 0.;

  /* allocate memory */
  e = gsl_vector_alloc(quad_num);
  r_last = gsl_vector_alloc(r_size);

  gsl_vector_memcpy(e,ortho_basis);
  normalise(wQuad, e);

  /* if not done, R matrix fills up below diagonal with .5 instead of 0 */
  gsl_vector_set_zero(ru);  

  while( !flag ){
    gsl_vector_memcpy(ortho_basis, e);
    gsl_vector_set_zero(r_last);

    modified_gs(r_last, ortho_basis, RB_space, wQuad, dim_RB);

    gsl_vector_add(ru, r_last);
    nrm_current = gsl_vector_get(r_last, dim_RB);

    gsl_vector_scale(ortho_basis, nrm_current);

    if( nrm_current/nrm_prev <= ortho_condition ) {
      nrm_prev = nrm_current;
      iter = iter + 1;
      gsl_vector_memcpy(e, ortho_basis);
    }
    else{ flag = 1; }

    nrm_current = normalisation(wQuad, ortho_basis);
    gsl_vector_set(ru, dim_RB, nrm_current);

    normalise(wQuad, ortho_basis);
  }

  gsl_vector_free(e);
  gsl_vector_free(r_last);
}

728 729 730 731 732 733
/* main functions */

/**
 * \brief Create a orthonormal basis set from a training set of real waveforms
 *
 * Given a \c gsl_matrix containing a training set of real waveforms (where the waveforms
734
 * are created at time or frequency steps seperated by \c delta) an orthonormal basis
735
 * will be generated using the greedy binning Algorithm 1 of \cite FGHKT2014 . The stopping
736
 * criteria for the algorithm is controlled by the \c tolerance value, which defined the
737
 * maximum residual between the current basis set (at a given iteration) and the training
738 739 740 741 742
 * set (and example tolerance is \f$10^{-12}\f$). In this function the training set will be
 * normalised, so the input \c TS will be modified.
 *
 * If the \c RBin value is \c NULL then a new reduced basis will be formed from the given
 * training set. However, if \c RBin already contains a previously produced basis, then this
743 744 745
 * basis will be enriched with bases if possible using the new training set.  <b>NOTE</b>: when
 * using  small tolerances enriching the basis in this way can lead to numerical precision issues,
 * so in general you should use \c LALInferenceEnrichREAL8Basis for enrichment.
746
 *
747 748 749
 * @param[out] RBin A \c REAL8Array to return the reduced basis.
 * @param[in] delta The time/frequency step(s) in the training set used to normalise the models.
 * This can be a vector containing just one value.
750
 * @param[in] tolerance The tolerance used as a stopping criteria for the basis generation.
751 752 753 754 755
 * @param[in] TS A \c REAL8Array matrix containing the complex training set, where the number of
 * waveforms in the training set is given by the rows and the waveform points by the columns. This
 * will be modified by this function, as the training set will be normalised.
 * @param[out] greedypoints A \c UINT4Vector to return the indices of the training set rows that
 * have been used to form the reduced basis.
756
 *
757
 * @return A \c REAL8 with the maximum projection error for the final reduced basis.
758 759
 *
 * \sa LALInferenceEnrichREAL8Basis
760
 */
761
REAL8 LALInferenceGenerateREAL8OrthonormalBasis(REAL8Array **RBin,
762
                                                const REAL8Vector *delta,
763
                                                REAL8 tolerance,
764 765 766 767
                                                REAL8Array **TS,
                                                UINT4Vector **greedypoints){
  REAL8Array *ts = NULL;
  ts = *TS; // pointer to training set
768

769
  size_t cols = ts->dimLength->data[1], rows = ts->dimLength->data[0];
770

771
  /* view of training set */
772
  gsl_matrix_view TSview;
773
  XLAL_CALLGSL( TSview = gsl_matrix_view_array((double *)ts->data, rows, cols) );
774

775
  size_t max_RB = rows;
776

777 778 779 780
  //REAL8 greedy_err[max_RB];           /* approximate error */
  UINT4Vector *gpts = NULL;
  gpts = XLALCreateUINT4Vector(max_RB); /* selected greedy points (row selection) */
  *greedypoints = gpts;
781

782 783 784
  REAL8 worst_err;          /* errors in greedy sweep */
  UINT4 worst_app = 0;      /* worst error stored */
  REAL8 tmpc;               /* worst error temp */
785

786 787 788 789 790
  gsl_vector *ts_el, *last_rb, *ortho_basis, *ru;
  gsl_matrix *R_matrix;
  REAL8 A_row_norms2[rows];              // || A(i,:) ||^2
  REAL8 projection_norms2[rows];
  REAL8 errors[rows];                    // approximation errors at i^{th} sweep
791

792 793 794
  REAL8Array *RB = NULL;
  UINT4Vector *dims = NULL;
  gsl_matrix_view RBview;
795

796 797 798 799 800
  /* this memory should be freed here */
  ts_el         = gsl_vector_alloc(cols);
  last_rb       = gsl_vector_alloc(cols);
  ortho_basis   = gsl_vector_alloc(cols);
  ru            = gsl_vector_alloc(max_RB);
801

802 803 804
  //project_coeff = gsl_matrix_alloc(max_RB,rows);
  REAL8 projection_coeff;
  R_matrix = gsl_matrix_alloc(max_RB, max_RB);
805

806 807
  /* initialise projection norms with zeros */
  for(size_t i=0; i<rows; ++i){ projection_norms2[i] = 0; }
808

809 810
  gsl_vector_view deltaview;
  XLAL_CALLGSL( deltaview = gsl_vector_view_array(delta->data, delta->length) );
811

812 813
  /* normalise the training set */
  normalise_training_set(&deltaview.vector, &TSview.matrix);
814

815 816 817 818 819
  /* compute norm of each training space element */
  for(size_t i=0; i<rows; ++i) {
    gsl_matrix_get_row(ts_el, &TSview.matrix, i);
    A_row_norms2[i] = normalisation(&deltaview.vector, ts_el);
  }
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
  /* initialize algorithm with first training set value */
  dims = XLALCreateUINT4Vector( 2 );
  dims->data[0] = 1; /* one row */
  dims->data[1] = cols;
  RB = XLALCreateREAL8Array( dims );
  *RBin = RB;

  XLAL_CALLGSL( RBview = gsl_matrix_view_array((double*)RB->data, 1, cols) );
  gsl_matrix_get_row(ts_el, &TSview.matrix, 0);
  gsl_matrix_set_row(&RBview.matrix, 0, ts_el);
  tmpc = 1.;
  gsl_matrix_set(R_matrix, 0, 0, tmpc);  // assumes normalized solutions

  gpts->data[0] = 0;
  UINT4 dim_RB          = 1;
  //greedy_err[0]    = 1.0;

  /* loop to find reduced basis */
  while( 1 ){
    gsl_matrix_get_row(last_rb, &RBview.matrix, dim_RB-1); /* previous basis */

    /* Compute overlaps of pieces of training set with rb_new */
    for(size_t i = 0; i < rows; i++){
      gsl_matrix_get_row(ts_el, &TSview.matrix, i);
      projection_coeff = weighted_dot_product(&deltaview.vector, last_rb, ts_el);
      projection_norms2[i] += (projection_coeff*projection_coeff);
      errors[i] = A_row_norms2[i] - projection_norms2[i];
848 849
    }

850 851 852 853 854 855
    /* find worst represented training set element, and add to basis */
    worst_err = 0.0;
    for(size_t i = 0; i < rows; i++) {
      if(worst_err < errors[i]) {
        worst_err = errors[i];
        worst_app = i;
856
      }
857 858
    }

859 860 861 862 863 864 865 866 867 868 869 870
    // This block exists for cases when the reduced basis is being used during enrichment.
    // It exists to make sure that all new training elements get added to the reduced basis
    // even if their errors are very small.
    if ( tolerance == 0. ){
      // check if worst_app is already in gpts
      UINT4 idxexists = 0;
      for(size_t i = 0; i < dim_RB; i++) {
        if ( worst_app == gpts->data[i] ){
          idxexists = 1;
          break;
        }
      }
871

872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890
      // if it is, then just find the first index that is not in gpts already
      if ( idxexists ){
        UINT4 newidxexists = 0;
        for (size_t i = 0; i < rows; i++){
          newidxexists = 0;
          for (size_t j = 0; j < dim_RB; j++){
            if ( i == gpts->data[j] ){
              newidxexists = 1;
              break;
            }
          }
          if ( !newidxexists ){
            worst_app = i;
            break;
          }
        }
      }
    }
    //////////// end check ////////////
891

892
    gpts->data[dim_RB] = worst_app;
893

894 895
    /* add worst approximated solution to basis set */
    gsl_matrix_get_row(ortho_basis, &TSview.matrix, worst_app);
896 897 898 899 900 901
    iterated_modified_gm(ru, ortho_basis, &RBview.matrix, &deltaview.vector, dim_RB); /* use IMGS */

    /* check normalisation of generated orthogonal basis is not NaN (cause by a new orthogonal basis
      having zero residual with the current basis) - if this is the case do not add the new basis. */
    if ( gsl_isnan(gsl_vector_get(ru, dim_RB)) ){ break; }

902 903
    /* add to reduced basis */
    dims->data[0] = dim_RB+1; /* add row */
904
    RB = XLALResizeREAL8Array( RB, dims );
905 906

    /* add on next basis */
907 908 909 910
    XLAL_CALLGSL( RBview = gsl_matrix_view_array((double*)RB->data, dim_RB+1, cols) );

    gsl_matrix_set_row(&RBview.matrix, dim_RB, ortho_basis);
    gsl_matrix_set_row(R_matrix, dim_RB, ru);
911

912 913 914 915
    ++dim_RB;

    /* decide if another greedy sweep is needed */
    if( (dim_RB == max_RB) || (worst_err < tolerance) || (rows == dim_RB) ){ break; }
916 917
  }

918 919 920 921 922 923 924 925
  gpts = XLALResizeUINT4Vector( gpts, dim_RB );

  XLALDestroyUINT4Vector(dims);
  gsl_vector_free(ts_el);
  gsl_vector_free(last_rb);
  gsl_vector_free(ortho_basis);
  gsl_vector_free(ru);
  gsl_matrix_free(R_matrix);
926

927
  return worst_err;
928 929 930 931 932 933
}


/**
 * \brief Create a orthonormal basis set from a training set of complex waveforms
 *
934 935
 * Given a \c gsl_matrix containing a training set \c TS of complex waveforms (where the waveforms
 * are created at time or frequency steps seperated by \c delta) an orthonormal basis
936
 * will be generated using the greedy binning Algorithm 1 of \cite FGHKT2014 . The stopping
937
 * criteria for the algorithm is controlled by the \c tolerance value, which is defined the
938
 * maximum residual between the current basis set (at a given iteration) and the training
939 940 941 942 943
 * set (and example tolerance is \f$10^{-12}\f$). In this function the training set will be
 * normalised, so the input \c TS will be modified.
 *
 * If the \c RBin value is \c NULL then a new reduced basis will be formed from the given
 * training set. However, if \c RBin already contains a previously produced basis, then this
944 945 946
 * basis will be enriched with bases if possible using the new training set. <b>NOTE</b>: when
 * using  small tolerances enriching the basis in this way can lead to numerical precision issues,
 * so in general you should use \c LALInferenceEnrichCOMPLEX16Basis for enrichment.
947 948 949 950 951 952
 *
 * Note that in this function we have to cast the \c COMPLEX16 array as a double to use
 * \c gsl_matrix_view_array, which assume that the data is passed as a double array with
 * memory laid out so that adjacent double memory blocks hold the corresponding real and
 * imaginary parts.
 *
953
 * @param[out] RBin A \c COMPLEX16Array to return the reduced basis.
954 955 956
 * @param[in] delta The time/frequency step(s) in the training set used to normalise the models.
 * This can be a vector containing just one value.
 * @param[in] tolerance The tolerance used as a stopping criteria for the basis generation.
957 958 959 960 961
 * @param[in] TS A \c COMPLEX16Array matrix containing the complex training set, where the number
 * of waveforms in the training set is given by the rows and the waveform points by the columns.
 * This will be modified by this function, as the training set will be normalised.
 * @param[out] greedypoints A \c UINT4Vector to return the indices of the training set rows that
 * have been used to form the reduced basis.
962
 *
963
 * @return A \c REAL8 with the maximum projection error for the final reduced basis.
964 965
 *
 * \sa LALInferenceEnrichCOMPLEX16Basis
966
 */
967
REAL8 LALInferenceGenerateCOMPLEX16OrthonormalBasis(COMPLEX16Array **RBin,
968
                                                    const REAL8Vector *delta,
969
                                                    REAL8 tolerance,
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
                                                    COMPLEX16Array **TS,
                                                    UINT4Vector **greedypoints){
  COMPLEX16Array *ts = NULL;
  ts = *TS; // pointer to training set

  size_t cols = ts->dimLength->data[1], rows = ts->dimLength->data[0];

  /* view of training set */
  gsl_matrix_complex_view TSview;
  XLAL_CALLGSL( TSview = gsl_matrix_complex_view_array((double *)ts->data, rows, cols) );

  size_t max_RB = rows;

  UINT4Vector *gpts = NULL;
  gpts = XLALCreateUINT4Vector(max_RB); /* selected greedy points (row selection) */
  *greedypoints = gpts;

  REAL8 worst_err;          /* errors in greedy sweep */
  UINT4 worst_app = 0;      /* worst error stored */
  gsl_complex tmpc;         /* worst error temp */

  gsl_vector_complex *ts_el, *last_rb, *ortho_basis, *ru;
  gsl_matrix_complex *R_matrix;
  REAL8 A_row_norms2[rows];              // || A(i,:) ||^2
  REAL8 projection_norms2[rows];
  REAL8 errors[rows];                    // approximation errors at i^{th} sweep
996 997 998

  COMPLEX16Array *RB = NULL;
  UINT4Vector *dims = NULL;
999
  gsl_matrix_complex_view RBview;
1000 1001 1002 1003 1004 1005
  
  /* this memory should be freed here */
  ts_el         = gsl_vector_complex_alloc(cols);
  last_rb       = gsl_vector_complex_alloc(cols);
  ortho_basis   = gsl_vector_complex_alloc(cols);
  ru            = gsl_vector_complex_alloc(max_RB);
1006

1007 1008 1009
  //project_coeff = gsl_matrix_complex_alloc(max_RB,rows);
  gsl_complex projection_coeff;
  R_matrix = gsl_matrix_complex_alloc(max_RB, max_RB);
1010

1011 1012
  /* initialise projection norms with zeros */
  for(size_t i=0; i<rows; ++i){ projection_norms2[i] = 0; }
1013

1014 1015
  gsl_vector_view deltaview;
  XLAL_CALLGSL( deltaview = gsl_vector_view_array(delta->data, delta->length) );
1016 1017
  
  /* normalise the training set */
1018 1019
  complex_normalise_training_set(&deltaview.vector, &TSview.matrix);

1020 1021 1022 1023 1024
  /* compute norm of each training space element */
  for(size_t i=0; i<rows; ++i) {
    gsl_matrix_complex_get_row(ts_el, &TSview.matrix, i);
    A_row_norms2[i] = complex_normalisation(&deltaview.vector, ts_el);
  }
1025

1026
  /* initialize algorithm with first training set value */
1027
  dims = XLALCreateUINT4Vector( 2 );
1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052
  dims->data[0] = 1; /* one row */
  dims->data[1] = cols;
  RB = XLALCreateCOMPLEX16Array( dims );
  *RBin = RB;

  XLAL_CALLGSL( RBview = gsl_matrix_complex_view_array((double*)RB->data, 1, cols) );
  gsl_matrix_complex_get_row(ts_el, &TSview.matrix, 0);
  gsl_matrix_complex_set_row(&RBview.matrix, 0, ts_el);
  GSL_SET_COMPLEX(&tmpc, 1.0, 0.0);
  gsl_matrix_complex_set(R_matrix, 0, 0, tmpc);  // assumes normalized solutions

  gpts->data[0] = 0;
  UINT4 dim_RB          = 1;

  /* loop to find reduced basis */
  while( 1 ){
    gsl_matrix_complex_get_row(last_rb, &RBview.matrix, dim_RB-1); /* previous basis */

    /* Compute overlaps of pieces of training set with rb_new */
    for(size_t i = 0; i < rows; i++){
      gsl_matrix_complex_get_row(ts_el, &TSview.matrix, i);
      projection_coeff = complex_weighted_dot_product(&deltaview.vector, last_rb, ts_el);
      projection_norms2[i] += (projection_coeff.dat[0]*projection_coeff.dat[0] + projection_coeff.dat[1]*projection_coeff.dat[1]);
      errors[i] = A_row_norms2[i] - projection_norms2[i];
    }
1053

1054 1055 1056 1057 1058 1059 1060 1061
    /* find worst represented training set element, and add to basis */
    worst_err = 0.0;
    for(size_t i = 0; i < rows; i++) {
      if(worst_err < errors[i]) {
        worst_err = errors[i];
        worst_app = i;
      }
    }
1062

1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074
    // This block exists for cases when the reduced basis is being used during enrichment.
    // It exists to make sure that all new training elements get added to the reduced basis
    // even if their errors are very small.
    if ( tolerance == 0. ){
      // check if worst_app is already in gpts
      UINT4 idxexists = 0;
      for(size_t i = 0; i < dim_RB; i++) {
        if ( worst_app == gpts->data[i] ){
          idxexists = 1;
          break;
        }
      }
1075

1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094
      // if it is, then just find the first index that is not in gpts already
      if ( idxexists ){
        UINT4 newidxexists = 0;
        for (size_t i = 0; i < rows; i++){
          newidxexists = 0;
          for (size_t j = 0; j < dim_RB; j++){
            if ( i == gpts->data[j] ){
              newidxexists = 1;
              break;
            }
          }
          if ( !newidxexists ){
            worst_app = i;
            break;
          }
        }
      }
    }
    //////////// end check ////////////
1095

1096
    gpts->data[dim_RB] = worst_app;
1097

1098 1099 1100
    /* add worst approximated solution to basis set */
    gsl_matrix_complex_get_row(ortho_basis, &TSview.matrix, worst_app);
    iterated_modified_gm_complex(ru, ortho_basis, &RBview.matrix, &deltaview.vector, dim_RB); /* use IMGS */
1101 1102 1103 1104 1105

    /* check normalisation of generated orthogonal basis is not NaN (cause by a new orthogonal basis
      having zero residual with the current basis) - if this is the case do not add the new basis. */
    if ( gsl_isnan(GSL_REAL(gsl_vector_complex_get(ru, dim_RB))) ){ break; }

1106 1107 1108
    /* add to reduced basis */
    dims->data[0] = dim_RB+1; /* add row */
    RB = XLALResizeCOMPLEX16Array( RB, dims );
1109

1110 1111
    /* add on next basis */
    XLAL_CALLGSL( RBview = gsl_matrix_complex_view_array((double*)RB->data, dim_RB+1, cols) );
1112

1113 1114
    gsl_matrix_complex_set_row(&RBview.matrix, dim_RB, ortho_basis);
    gsl_matrix_complex_set_row(R_matrix, dim_RB, ru);
1115

1116
    ++dim_RB;
1117

1118 1119 1120
    /* decide if another greedy sweep is needed */
    if( (dim_RB == max_RB) || (worst_err < tolerance) || (rows == dim_RB) ){ break; }
  }
1121

1122
  gpts = XLALResizeUINT4Vector( gpts, dim_RB );
1123

1124 1125 1126 1127 1128 1129
  XLALDestroyUINT4Vector(dims);
  gsl_vector_complex_free(ts_el);
  gsl_vector_complex_free(last_rb);
  gsl_vector_complex_free(ortho_basis);
  gsl_vector_complex_free(ru);
  gsl_matrix_complex_free(R_matrix);
1130

1131 1132
  return worst_err;
}
1133

1134

1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163
/**
 * \brief Validate the real reduced basis against another set of waveforms
 *
 * This function projects a set of waveforms onto the reduced basis and
 * checks that the residuals are within a given tolerance. It returns
 * the projection errors.
 *
 * Note that the projection error returned are the square of the residual
 * errors, as is used as the criterion for adding new bases in the reduced
 * basis generation function \c LALInferenceGenerateREAL8OrthonormalBasis.
 * This is different to the \c validation.cpp code in <a href="https://bitbucket.org/sfield83/greedycpp">greedycpp</a>,
 * which returns the square root of the value we return.
 *
 * @param[out] projerr The projection errors for each test waveform.
 * @param[in] delta The time/frequency step(s) in the training set used to normalise the models.
 * This can be a vector containing just one value.
 * @param[in] RB The reduced basis set
 * @param[in] testmodels The set of waveform models to project onto the basis (these will
 * be changed by this function, as the waveforms will get normalised).
 *
 * \sa LALInferenceTestREAL8OrthonormalBasis
 */
void LALInferenceValidateREAL8OrthonormalBasis(REAL8Vector **projerr,
                                               const REAL8Vector *delta,
                                               const REAL8Array *RB,
                                               REAL8Array **testmodels){
  XLAL_CHECK_VOID( delta != NULL, XLAL_EFUNC, "Vector of 'delta' values is NULL!" );
  XLAL_CHECK_VOID( RB != NULL, XLAL_EFUNC, "Reduced basis set array is NULL!" );
  XLAL_CHECK_VOID( RB->dimLength->length == 2, XLAL_EFUNC, "Reduced basis set array must have only two dimensions" );
1164

1165 1166
  REAL8Array *tm = NULL;
  tm = *testmodels;
1167

1168 1169 1170
  size_t dlength = RB->dimLength->data[1], nts = tm->dimLength->data[0];
  size_t k = 0;
  gsl_vector *r_tmp, *testrow;
1171

1172 1173 1174 1175 1176 1177
  /* normalise the test set */
  gsl_vector_view deltaview;
  XLAL_CALLGSL( deltaview = gsl_vector_view_array(delta->data, delta->length) );
  gsl_matrix_view testmodelsview;
  XLAL_CALLGSL( testmodelsview = gsl_matrix_view_array(tm->data, nts, dlength) );
  normalise_training_set(&deltaview.vector, &testmodelsview.matrix);
1178

1179 1180
  gsl_matrix_view RBview;
  RBview = gsl_matrix_view_array(RB->data, RB->dimLength->data[0], dlength);
1181

1182 1183
  XLAL_CALLGSL( testrow = gsl_vector_alloc(dlength) );
  XLAL_CALLGSL( r_tmp = gsl_vector_alloc(RB->dimLength->data[0]) );
1184

1185 1186 1187
  REAL8Vector *pe = NULL;
  pe = XLALCreateREAL8Vector( nts );
  *projerr = pe;
1188

1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209
  for ( k = 0; k < nts; k++ ){
    pe->data[k] = 0.;
    REAL8 r_tmp_nrm = 0.;

    XLAL_CALLGSL( gsl_matrix_get_row(testrow, &testmodelsview.matrix, k) );

    REAL8 nrm = normalisation(&deltaview.vector, testrow); // normalisation (should be 1 as test models are normalised)

    // un-scale testrow
    if ( delta->length == 1 ){
      XLAL_CALLGSL( gsl_vector_scale(testrow, delta->data[0]) );
    }
    else if ( testrow->size == deltaview.vector.size ){
      XLAL_CALLGSL( gsl_vector_mul(testrow, &deltaview.vector) );
    }
    else{
      XLAL_ERROR_VOID(XLAL_EFUNC, "Vector of weights must either contain a single value, or be the same length as the other input vectors.");
    }

    // get projections
    XLAL_CALLGSL( gsl_blas_dgemv( CblasNoTrans, 1., &RBview.matrix, testrow, 0., r_tmp ) );
1210

1211 1212 1213 1214 1215 1216 1217 1218
    r_tmp_nrm = gsl_blas_dnrm2( r_tmp );
    pe->data[k] = nrm - r_tmp_nrm*r_tmp_nrm;

    if ( pe->data[k] < 0. ) { pe->data[k] = 1.0e-16; } // floating point error can trigger this
  }

  XLAL_CALLGSL( gsl_vector_free( testrow ) );
  XLAL_CALLGSL( gsl_vector_free( r_tmp ) );
1219 1220 1221 1222
}


/**
1223
 * \brief Test the reduced basis against another set of waveforms
1224 1225 1226 1227 1228 1229
 *
 * This function projects a set of waveforms onto the reduced basis and
 * checks that the residuals are within a given tolerance
 *
 * @param[in] delta The time/frequency step(s) in the training set used to normalise the models.
 * This can be a vector containing just one value.
1230
 * @param[in] tolerance The allowed (squared) residual tolerence for the test waveforms
1231 1232 1233
 * @param[in] RB The reduced basis set
 * @param[in] testmodels The set of waveform models to project onto the basis
 *
1234
 * @return Returns \c XLAL_SUCCESS if all test waveforms meet the tolerance
1235 1236
 *
 * \sa LALInferenceValidateREAL8OrthonormalBasis
1237
 */
1238
INT4 LALInferenceTestREAL8OrthonormalBasis(const REAL8Vector *delta,
1239
                                           REAL8 tolerance,
1240 1241
                                           const REAL8Array *RB,
                                           REAL8Array **testmodels){
1242 1243 1244 1245 1246
  XLAL_CHECK( delta != NULL, XLAL_EFUNC, "Vector of 'delta' values is NULL!" );
  XLAL_CHECK( RB != NULL, XLAL_EFUNC, "Reduced basis set array is NULL!" );
  XLAL_CHECK( RB->dimLength->length == 2, XLAL_EFUNC, "Reduced basis set array must have only two dimensions" );
  XLAL_CHECK( tolerance > 0, XLAL_EFUNC, "Tolerance is less than, or equal to, zero!" );

1247 1248
  REAL8Array *tm = NULL;
  tm = *testmodels;
1249

1250
  size_t nts = tm->dimLength->data[0], k = 0;
1251

1252
  REAL8Vector *projerr = NULL;
1253

1254
  LALInferenceValidateREAL8OrthonormalBasis(&projerr, delta, RB, testmodels);
1255

1256
  for ( k = 0; k < nts; k++ ){
1257
    /* check projection error against tolerance */
1258 1259 1260 1261
    if ( projerr->data[k] > tolerance ) {
      XLALDestroyREAL8Vector( projerr );
      return XLAL_FAILURE;
    }
1262 1263
  }

1264 1265
  XLALDestroyREAL8Vector( projerr );

1266 1267 1268
  return XLAL_SUCCESS;
}

1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279
/**
 * \brief Expand the real training waveforms with ones from a set of new training waveforms
 *
 * Expands an original set of training waveforms with waveforms from a new set a training waveforms
 * if, when projected onto the reduced basis, those new waveforms are greater than the given tolerance.
 * The reduced basis will then be recomputed using the expanded training set.
 *
 * @param[in] delta The time/frequency step(s) in the training set used to normalise the models.
 * This can be a vector containing just one value.
 * @param[in] tolerance The allowed residual tolerence for the test waveforms
 * @param[in] RB The reduced basis set
1280 1281 1282 1283 1284 1285
 * @param[in] greedypoints The rows in \c testmodels that formed the reduced basis
 * @param[in] testmodels The set of waveform models used to produce the reduced basis (already normalised)
 * @param[in] testmodelsnew A new set of waveforms to project onto the current reduced basis (these will
 * be changed by this function, as the waveforms will get normalised when pass to
 * \c LALInferenceValidateREAL8OrthonormalBasis, and editted to contain the new set of training waveforms
 * from which the enriched basis was formed).
1286 1287 1288
 *
 * @return Returns the maximum projection error of the new basis set
 */
1289
REAL8 LALInferenceEnrichREAL8Basis(const REAL8Vector *delta,
1290 1291
                                   REAL8 tolerance,
                                   REAL8Array **RB,
1292 1293 1294
                                   UINT4Vector **greedypoints,
                                   const REAL8Array *testmodels,
                                   REAL8Array **testmodelsnew){
1295 1296 1297
  XLAL_CHECK( delta != NULL, XLAL_EFUNC, "Vector of 'delta' values is NULL!" );
  XLAL_CHECK( tolerance > 0, XLAL_EFUNC, "Tolerance is less than, or equal to, zero!" );

1298 1299
  size_t dlength = testmodels->dimLength->data[1]; // length of training vectors
  size_t k = 0;
1300 1301 1302 1303 1304

  REAL8 maxprojerr = 0.;

  REAL8Array *RBin = NULL;
  RBin = *RB;
1305 1306
  UINT4Vector *gp = NULL;
  gp = *greedypoints;
1307
  REAL8Array *tm = NULL;
1308
  tm = *testmodelsnew;
1309

1310
  size_t nts = tm->dimLength->data[0]; // number of new training vectors
1311 1312

  XLAL_CHECK_REAL8( dlength == tm->dimLength->data[1], XLAL_EFUNC, "New training set contains waveforms of different length to the original set!" );
1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348
  
  gsl_vector_view deltaview;
  XLAL_CALLGSL( deltaview = gsl_vector_view_array(delta->data, delta->length) );
  
  // get projection errors of new test model set against the reduced basis
  REAL8Vector *projerr = NULL;
  LALInferenceValidateREAL8OrthonormalBasis( &projerr, delta, *RB, testmodelsnew );

  // check for errors outside the tolerance value, and reduced new training set accordingly
  size_t keepcounter = 0;
  gsl_matrix_view tmview;
  tmview = gsl_matrix_view_array(tm->data, nts, dlength);
  for ( k = 0; k < projerr->length; k++ ){
    if ( projerr->data[k] > tolerance ){
      if ( keepcounter != k ){
        gsl_vector_view row;
        row = gsl_matrix_row( &tmview.matrix, k );

        // un-scale row (so that it works with the reduced basis generation function below)
        if ( delta->length == 1 ){
          XLAL_CALLGSL( gsl_vector_scale(&row.vector, delta->data[0]) );
        }
        else if ( row.vector.size == delta->length ){
          XLAL_CALLGSL( gsl_vector_mul(&row.vector, &deltaview.vector) );
        }
        else{
          XLAL_ERROR_REAL8(XLAL_EFUNC, "Vector of weights must either contain a single value, or be the same length as the other input vectors.");
        }
        
        gsl_matrix_set_row( &tmview.matrix, keepcounter, &row.vector ); // copy row to earlier in the matrix
      }
      keepcounter++;
    }
  }

  XLALDestroyREAL8Vector( projerr );
1349

1350 1351
  if ( keepcounter == 0 ){ return maxprojerr; } // nothing was above the tolerance

1352
  // add vectors used to produce the original reduced basis
1353
  UINT4Vector *dims = XLALCreateUINT4Vector( 2 );
1354
  dims->data[0] = RBin->dimLength->data[0] + keepcounter;
1355
  dims->data[1] = dlength;
1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367
  tm = XLALResizeREAL8Array( tm, dims );
  tmview = gsl_matrix_view_array(tm->data, dims->data[0], dlength);

  gsl_matrix_view otmview;
  otmview = gsl_matrix_view_array(testmodels->data, testmodels->dimLength->data[0], dlength); // view of original training set
  for ( k = 0; k < RBin->dimLength->data[0]; k++ ){
    gsl_vector_view row;
    row = gsl_matrix_row( &otmview.matrix, gp->data[k] );

    // un-scale row (so that it works with the reduced basis generation function below)
    if ( delta->length == 1 ){
      XLAL_CALLGSL( gsl_vector_scale(&row.vector, delta->data[0]) );
1368
    }
1369 1370 1371 1372 1373
    else if ( row.vector.size == delta->length ){
      XLAL_CALLGSL( gsl_vector_mul(&row.vector, &deltaview.vector) );
    }
    else{
      XLAL_ERROR_REAL8(XLAL_EFUNC, "Vector of weights must either contain a single value, or be the same length as the other input vectors.");
1374 1375
    }

1376
    gsl_matrix_set_row( &tmview.matrix, keepcounter + k, &row.vector );
1377 1378 1379
  }

  /* recalculate the reduced basis using the updated training set */
1380 1381 1382 1383 1384 1385 1386
  XLALDestroyREAL8Array( *RB ); /* remove current basis, so it is rebuilt using the entire new training set */
  *RB = NULL;
  XLALDestroyUINT4Vector( *greedypoints );
  *greedypoints = NULL;
  maxprojerr = LALInferenceGenerateREAL8OrthonormalBasis(RB, delta, 0.0, &tm, greedypoints); // set tolerance to 0, so that points are always added

  XLALDestroyUINT4Vector( dims );
1387 1388 1389 1390 1391

  return maxprojerr;
}


1392
/**
1393
 * \brief Validate the complex reduced basis against another set of waveforms
1394 1395
 *
 * This function projects a set of waveforms onto the reduced basis and
1396 1397 1398 1399 1400 1401 1402 1403
 * checks that the residuals are within a given tolerance. It returns
 * the projection errors.
 *
 * Note that the projection error returned are the square of the residual
 * errors, as is used as the criterion for adding new bases in the reduced
 * basis generation function \c LALInferenceGenerateCOMPLEX16OrthonormalBasis.
 * This is different to the \c validation.cpp code in <a href="https://bitbucket.org/sfield83/greedycpp">greedycpp</a>,
 * which returns the square root of the value we return.
1404
 *
1405
 * @param[out] projerr The projection errors (square of the residual) for each test waveform.
1406 1407 1408
 * @param[in] delta The time/frequency step(s) in the training set used to normalise the models.
 * This can be a vector containing just one value.
 * @param[in] RB The reduced basis set
1409 1410
 * @param[in] testmodels The set of waveform models to project onto the basis (these will
 * be changed by this function, as the waveforms will get normalised).
1411
 *
1412
 * \sa LALInferenceTestCOMPLEX16OrthonormalBasis
1413
 */
1414 1415 1416 1417 1418 1419 1420
void LALInferenceValidateCOMPLEX16OrthonormalBasis(REAL8Vector **projerr,
                                                   const REAL8Vector *delta,
                                                   const COMPLEX16Array *RB,
                                                   COMPLEX16Array **testmodels){
  XLAL_CHECK_VOID( delta != NULL, XLAL_EFUNC, "Vector of 'delta' values is NULL!" );
  XLAL_CHECK_VOID( RB != NULL, XLAL_EFUNC, "Reduced basis is NULL!" );
  XLAL_CHECK_VOID( RB->dimLength->length == 2, XLAL_EFUNC, "Reduced basis set does not have 2 dimensions!" );
1421

1422 1423 1424 1425 1426 1427
  COMPLEX16Array *tm = NULL;
  tm = *testmodels;

  size_t dlength = RB->dimLength->data[1], nts = tm->dimLength->data[0];
  size_t k = 0, j = 0;
  gsl_vector_complex *r_tmp, *testrow;
1428 1429

  /* normalise the test set */
1430 1431 1432
  gsl_vector_view deltaview;
  XLAL_CALLGSL( deltaview = gsl_vector_view_array(delta->data, delta->length) );
  gsl_matrix_complex_view testmodelsview;
1433
  XLAL_CALLGSL( testmodelsview = gsl_matrix_complex_view_array((double *)tm->data, nts, dlength) );
1434 1435 1436 1437
  complex_normalise_training_set(&deltaview.vector, &testmodelsview.matrix);

  gsl_matrix_complex_view RBview;
  RBview = gsl_matrix_complex_view_array((double *)RB->data, RB->dimLength->data[0], dlength);
1438

1439 1440 1441 1442 1443 1444 1445
  XLAL_CALLGSL( testrow = gsl_vector_complex_alloc(dlength) );
  XLAL_CALLGSL( r_tmp = gsl_vector_complex_alloc(RB->dimLength->data[0]) );

  REAL8Vector *pe = NULL;
  pe = XLALCreateREAL8Vector( nts );
  *projerr = pe;

1446 1447
  /* get projection errors for each test model */
  for ( k = 0; k < nts; k++ ){
1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467
    pe->data[k] = 0.;
    REAL8 r_tmp_nrm = 0.;

    XLAL_CALLGSL( gsl_matrix_complex_get_row(testrow, &testmodelsview.matrix, k) );

    REAL8 nrm = complex_normalisation(&deltaview.vector, testrow); // normalisation (should be 1 as test models are normalised)

    // un-scale testrow (gsl_vector_complex_mul doesn't work when delta is a real vector)
    if ( delta->length == 1 || testrow->size == deltaview.vector.size ){
      for ( j = 0; j < dlength; j++ ){
        gsl_complex ctmp;
        XLAL_CALLGSL( ctmp = gsl_vector_complex_get( testrow, j ) );
        if ( delta->length == 1 ){
          XLAL_CALLGSL( ctmp = gsl_complex_mul_real( ctmp, delta->data[0] ) );
        }
        else{
          XLAL_CALLGSL( ctmp = gsl_complex_mul_real( ctmp, delta->data[j] ) );
        }
        XLAL_CALLGSL( gsl_vector_complex_set( testrow, j, ctmp ) );
      }
1468
    }
1469 1470 1471 1472 1473 1474
    else{
      XLAL_ERROR_VOID( XLAL_EFUNC, "Vector of weights must either contain a single value, or be the same length as the other input vectors." );
    }
    
    // get complex conjugate
    for ( j = 0; j < testrow->size; j++ ){ testrow->data[2*j+1] = -testrow->data[2*j+1]; }
1475

1476 1477
    // get projections
    XLAL_CALLGSL( gsl_blas_zgemv( CblasNoTrans, GSL_COMPLEX_ONE, &RBview.matrix, testrow, GSL_COMPLEX_ZERO, r_tmp ) );
1478

1479 1480
    r_tmp_nrm = gsl_blas_dznrm2(r_tmp);
    pe->data[k] = nrm - r_tmp_nrm*r_tmp_nrm;
1481

1482 1483 1484 1485 1486 1487 1488 1489 1490 1491 1492 1493 1494 1495 1496 1497 1498 1499 1500 1501 1502 1503 1504 1505 1506 1507 1508 1509 1510 1511 1512 1513 1514 1515 1516
    if ( pe->data[k] < 0. ) { pe->data[k] = 1.0e-16; } // floating point error can trigger this
  }

  XLAL_CALLGSL( gsl_vector_complex_free( testrow ) );
  XLAL_CALLGSL( gsl_vector_complex_free( r_tmp ) );
}


/**
 * \brief Test the reduced basis against another set of waveforms
 *
 * This function projects a set of waveforms onto the reduced basis and
 * checks that the residuals are within a given tolerance
 *
 * @param[in] delta The time/frequency step(s) in the training set used to normalise the models.
 * This can be a vector containing just one value.
 * @param[in] tolerance The allowed residual (squared) tolerence for the test waveforms
 * @param[in] RB The reduced basis set
 * @param[in] testmodels The set of waveform models to project onto the basis
 *
 * @return Returns \c XLAL_SUCCESS if all test waveforms meet the tolerance
 *
 * \sa LALInferenceValidateCOMPLEX16OrthonormalBasis
 */
INT4 LALInferenceTestCOMPLEX16OrthonormalBasis(const REAL8Vector *delta,
                                               REAL8 tolerance,
                                               const COMPLEX16Array *RB,
                                               COMPLEX16Array **testmodels){
  XLAL_CHECK( delta != NULL, XLAL_EFUNC, "Vector of 'delta' values is NULL!" );
  XLAL_CHECK( RB != NULL, XLAL_EFUNC, "Reduced basis is NULL!" );
  XLAL_CHECK( RB->dimLength->length == 2, XLAL_EFUNC, "Reduced basis set does not have 2 dimensions!" );
  XLAL_CHECK( tolerance > 0, XLAL_EFUNC, "Tolerance is less than, or equal to, zero!" );

  COMPLEX16Array *tm = NULL;
  tm = *testmodels;
1517

1518 1519 1520 1521 1522 1523 1524 1525 1526 1527 1528 1529
  size_t nts = tm->dimLength->data[0], k = 0;

  REAL8Vector *projerr = NULL;

  LALInferenceValidateCOMPLEX16OrthonormalBasis(&projerr, delta, RB, testmodels);

  for ( k = 0; k < nts; k++ ){
    /* check projection error against tolerance */
    if ( projerr->data[k] > tolerance ) {
      XLALDestroyREAL8Vector( projerr );
      return XLAL_FAILURE;
    }
1530 1531
  }

1532 1533
  XLALDestroyREAL8Vector( projerr );

1534 1535 1536 1537
  return XLAL_SUCCESS;
}


1538 1539 1540 1541 1542 1543 1544 1545 1546 1547 1548
/**
 * \brief Expand the complex training waveforms with ones from a set of new training waveforms
 *
 * Expands an original set of training waveforms with waveforms from a new set a training waveforms
 * if, when projected onto the reduced basis, those new waveforms are greater than the given tolerance.
 * The reduced basis will then be recomputed using the expanded training set.
 *
 * @param[in] delta The time/frequency step(s) in the training set used to normalise the models.
 * This can be a vector containing just one value.
 * @param[in] tolerance The allowed residual tolerence for the test waveforms
 * @param[in] RB The reduced basis set
1549
 * @param[in] greedypoints The rows in \c testmodels that formed the reduced basis
1550
 * @param[in] testmodels The set of waveform models used to produce the reduced basis
1551 1552 1553 1554
 * @param[in] testmodelsnew A new set of waveforms to project onto the current reduced basis (these will
 * be changed by this function, as the waveforms will get normalised when pass to
 * \c LALInferenceValidateCOMPLEX16OrthonormalBasis, and editted to contain the new set of training waveforms
 * from which the enriched basis was formed)
1555 1556 1557
 *
 * @return Returns the maximum projection error of the new basis set
 */
1558
REAL8 LALInferenceEnrichCOMPLEX16Basis(const REAL8Vector *delta,
1559 1560
                                       REAL8 tolerance,
                                       COMPLEX16Array **RB,
1561 1562 1563
                                       UINT4Vector **greedypoints,
                                       const COMPLEX16Array *testmodels,
                                       COMPLEX16Array **testmodelsnew){
1564 1565 1566
  XLAL_CHECK( delta != NULL, XLAL_EFUNC, "Vector of 'delta' values is NULL!" );
  XLAL_CHECK( tolerance > 0, XLAL_EFUNC, "Tolerance is less than, or equal to, zero!" );

1567 1568
  size_t dlength = testmodels->dimLength->data[1]; // length of training vectors
  size_t k = 0, j = 0;