libtpcmodel
 All Data Structures Files Functions Variables Typedefs Enumerations Enumerator Defines
Functions
nnls.c File Reference
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include "include/nnls.h"
#include "include/lss_lib.h"
Include dependency graph for nnls.c:

Go to the source code of this file.

Functions

int nnls (double **a, int m, int n, double *b, double *x, double *rnorm, double *wp, double *zzp, int *indexp)
int nnlsWght (int N, int M, double **A, double *b, double *weight)
int nnlsWghtSquared (int N, int M, double **A, double *b, double *sweight)

Function Documentation

int nnls ( double **  a,
int  m,
int  n,
double *  b,
double *  x,
double *  rnorm,
double *  wp,
double *  zzp,
int *  indexp 
)

Algorithm NNLS (Non-negative least-squares)

Given an m by n matrix A, and an m-vector B, computes an n-vector X, that solves the least squares problem A * X = B , subject to X>=0

Instead of pointers for working space, NULL can be given to let this function to allocate and free the required memory.

Returns:
Function returns 0 if succesful, 1, if iteration count exceeded 3*N, or 2 in case of invalid problem dimensions or memory allocation error.
Parameters:
aOn entry, a[ 0... N ][ 0 ... M ] contains the M by N matrix A. On exit, a[][] contains the product matrix Q*A, where Q is an m by n orthogonal matrix generated implicitly by this function.
mMatrix dimension m
nMatrix dimension n
bOn entry, b[] must contain the m-vector B. On exit, b[] contains Q*B
xOn exit, x[] will contain the solution vector
rnormOn exit, rnorm contains the Euclidean norm of the residual vector. If NULL is given, no rnorm is calculated
wpAn n-array of working space, wp[]. On exit, wp[] will contain the dual solution vector. wp[i]=0.0 for all i in set p and wp[i]<=0.0 for all i in set z. Can be NULL, which causes this algorithm to allocate memory for it.
zzpAn m-array of working space, zz[]. Can be NULL, which causes this algorithm to allocate memory for it.
indexpAn n-array of working space, index[]. Can be NULL, which causes this algorithm to allocate memory for it.

Definition at line 63 of file nnls.c.

References _lss_g1(), _lss_h12(), _lss_print_matrix(), and NNLS_TEST.

Referenced by do_nnls(), and ldp().

Here is the call graph for this function:

int nnlsWght ( int  N,
int  M,
double **  A,
double *  b,
double *  weight 
)

Algorithm for weighting the problem that is given to nnls-algorithm. Square roots of weights are used because in nnls the difference w*A-w*b is squared.

Returns:
Algorithm returns zero if successful, 1 if arguments are inappropriate.
Parameters:
NNNLS dimension N (nr of parameters)
MNNLS dimension M (nr of samples)
ANNLS matrix A
bNNLS vector B
weightWeights for each sample (array of length M)

Definition at line 307 of file nnls.c.

References M.

int nnlsWghtSquared ( int  N,
int  M,
double **  A,
double *  b,
double *  sweight 
)

Algorithm for weighting the problem that is given to nnls-algorithm. Square roots of weights are used because in nnls the difference w*A-w*b is squared. Here user must give squared weights; this makes calculation faster, when this function needs to be called many times.

Returns:
Algorithm returns zero if successful, 1 if arguments are inappropriate.
Parameters:
NNNLS dimension N (nr of parameters)
MNNLS dimension M (nr of samples)
ANNLS matrix A
bNNLS vector B
sweightSquared weights for each sample (array of length M)

Definition at line 355 of file nnls.c.

References M.