1c4762a1bSJed Brown /* XH: todo add cs1f.F90 and asjust makefile */ 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown Include "petsctao.h" so that we can use TAO solvers. Note that this 4c4762a1bSJed Brown file automatically includes libraries such as: 5c4762a1bSJed Brown petsc.h - base PETSc routines petscvec.h - vectors 6a5b23f4aSJose E. Roman petscsys.h - system routines petscmat.h - matrices 7c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 8c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 9c4762a1bSJed Brown 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown 12c4762a1bSJed Brown #include <petsctao.h> 13c4762a1bSJed Brown 14c4762a1bSJed Brown /* 15c4762a1bSJed Brown Description: Compressive sensing test example 1. 16c4762a1bSJed Brown 0.5*||Ax-b||^2 + lambda*||D*x||_1 17c4762a1bSJed Brown Xiang Huang: Nov 19, 2018 18c4762a1bSJed Brown 19c4762a1bSJed Brown Reference: None 20c4762a1bSJed Brown */ 21c4762a1bSJed Brown 22c4762a1bSJed Brown static char help[] = "Finds the least-squares solution to the under constraint linear model Ax = b, with L1-norm regularizer. \n\ 23c4762a1bSJed Brown A is a M*N real matrix (M<N), x is sparse. \n\ 24c4762a1bSJed Brown We find the sparse solution by solving 0.5*||Ax-b||^2 + lambda*||D*x||_1, where lambda (by default 1e-4) is a user specified weight.\n\ 25c4762a1bSJed Brown D is the K*N transform matrix so that D*x is sparse. By default D is identity matrix, so that D*x = x.\n"; 26c4762a1bSJed Brown 27c4762a1bSJed Brown #define M 3 28c4762a1bSJed Brown #define N 5 29c4762a1bSJed Brown #define K 4 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* User-defined application context */ 32c4762a1bSJed Brown typedef struct { 33c4762a1bSJed Brown /* Working space. linear least square: f(x) = A*x - b */ 34c4762a1bSJed Brown PetscReal A[M][N]; /* array of coefficients */ 35c4762a1bSJed Brown PetscReal b[M]; /* array of observations */ 36c4762a1bSJed Brown PetscReal xGT[M]; /* array of ground truth object, which can be used to compare the reconstruction result */ 37c4762a1bSJed Brown PetscReal D[K][N]; /* array of coefficients for 0.5*||Ax-b||^2 + lambda*||D*x||_1 */ 38c4762a1bSJed Brown PetscReal J[M][N]; /* dense jacobian matrix array. For linear least square, J = A. For nonlinear least square, it is different from A */ 39c4762a1bSJed Brown PetscInt idm[M]; /* Matrix row, column indices for jacobian and dictionary */ 40c4762a1bSJed Brown PetscInt idn[N]; 41c4762a1bSJed Brown PetscInt idk[K]; 42c4762a1bSJed Brown } AppCtx; 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* User provided Routines */ 45c4762a1bSJed Brown PetscErrorCode InitializeUserData(AppCtx *); 46c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec); 47c4762a1bSJed Brown PetscErrorCode FormDictionaryMatrix(Mat, AppCtx *); 48c4762a1bSJed Brown PetscErrorCode EvaluateFunction(Tao, Vec, Vec, void *); 49c4762a1bSJed Brown PetscErrorCode EvaluateJacobian(Tao, Vec, Mat, Mat, void *); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /*--------------------------------------------------------------------*/ 52d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 53d71ae5a4SJacob Faibussowitsch { 54c4762a1bSJed Brown Vec x, f; /* solution, function f(x) = A*x-b */ 55c4762a1bSJed Brown Mat J, D; /* Jacobian matrix, Transform matrix */ 56c4762a1bSJed Brown Tao tao; /* Tao solver context */ 57c4762a1bSJed Brown PetscInt i; /* iteration information */ 58c4762a1bSJed Brown PetscReal hist[100], resid[100]; 59c4762a1bSJed Brown PetscInt lits[100]; 60c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 61c4762a1bSJed Brown 62327415f7SBarry Smith PetscFunctionBeginUser; 639566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* Allocate solution and vector function vectors */ 669566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &x)); 679566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, M, &f)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Allocate Jacobian and Dictionary matrix. */ 709566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, M, N, NULL, &J)); 719566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, K, N, NULL, &D)); /* XH: TODO: dense -> sparse/dense/shell etc, do it on fly */ 72c4762a1bSJed Brown 73c4762a1bSJed Brown for (i = 0; i < M; i++) user.idm[i] = i; 74c4762a1bSJed Brown for (i = 0; i < N; i++) user.idn[i] = i; 75c4762a1bSJed Brown for (i = 0; i < K; i++) user.idk[i] = i; 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 789566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF, &tao)); 799566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBRGN)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* User set application context: A, D matrice, and b vector. */ 829566063dSJacob Faibussowitsch PetscCall(InitializeUserData(&user)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Set initial guess */ 859566063dSJacob Faibussowitsch PetscCall(FormStartingPoint(x)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* Fill the content of matrix D from user application Context */ 889566063dSJacob Faibussowitsch PetscCall(FormDictionaryMatrix(D, &user)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Bind x to tao->solution. */ 919566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 92c4762a1bSJed Brown /* Bind D to tao->data->D */ 939566063dSJacob Faibussowitsch PetscCall(TaoBRGNSetDictionaryMatrix(tao, D)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* Set the function and Jacobian routines. */ 969566063dSJacob Faibussowitsch PetscCall(TaoSetResidualRoutine(tao, f, EvaluateFunction, (void *)&user)); 979566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianResidualRoutine(tao, J, J, EvaluateJacobian, (void *)&user)); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* Check for any TAO command line arguments */ 1009566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 101c4762a1bSJed Brown 1029566063dSJacob Faibussowitsch PetscCall(TaoSetConvergenceHistory(tao, hist, resid, 0, lits, 100, PETSC_TRUE)); 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* Perform the Solve */ 1059566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* XH: Debug: View the result, function and Jacobian. */ 1089566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "-------- result x, residual f=A*x-b, and Jacobian=A. -------- \n")); 1099566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_SELF)); 1109566063dSJacob Faibussowitsch PetscCall(VecView(f, PETSC_VIEWER_STDOUT_SELF)); 1119566063dSJacob Faibussowitsch PetscCall(MatView(J, PETSC_VIEWER_STDOUT_SELF)); 1129566063dSJacob Faibussowitsch PetscCall(MatView(D, PETSC_VIEWER_STDOUT_SELF)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Free TAO data structures */ 1159566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* Free PETSc data structures */ 1189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&f)); 1209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 122c4762a1bSJed Brown 1239566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 124b122ec5aSJacob Faibussowitsch return 0; 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown /*--------------------------------------------------------------------*/ 128d71ae5a4SJacob Faibussowitsch PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr) 129d71ae5a4SJacob Faibussowitsch { 130c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 131c4762a1bSJed Brown PetscInt m, n; 132c4762a1bSJed Brown const PetscReal *x; 133c4762a1bSJed Brown PetscReal *b = user->b, *f; 134c4762a1bSJed Brown 135c4762a1bSJed Brown PetscFunctionBegin; 1369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 1379566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 138c4762a1bSJed Brown 139a5b23f4aSJose E. Roman /* Even for linear least square, we do not direct use matrix operation f = A*x - b now, just for future modification and compatibility for nonlinear least square */ 140c4762a1bSJed Brown for (m = 0; m < M; m++) { 141c4762a1bSJed Brown f[m] = -b[m]; 142ad540459SPierre Jolivet for (n = 0; n < N; n++) f[m] += user->A[m][n] * x[n]; 143c4762a1bSJed Brown } 1449566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1459566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 146ca0c957dSBarry Smith PetscLogFlops(2.0 * M * N); 147c4762a1bSJed Brown PetscFunctionReturn(0); 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 150c4762a1bSJed Brown /*------------------------------------------------------------*/ 151c4762a1bSJed Brown /* J[m][n] = df[m]/dx[n] */ 152d71ae5a4SJacob Faibussowitsch PetscErrorCode EvaluateJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr) 153d71ae5a4SJacob Faibussowitsch { 154c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 155c4762a1bSJed Brown PetscInt m, n; 156c4762a1bSJed Brown const PetscReal *x; 157c4762a1bSJed Brown 158c4762a1bSJed Brown PetscFunctionBegin; 1599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); /* not used for linear least square, but keep for future nonlinear least square) */ 160c4762a1bSJed Brown /* XH: TODO: For linear least square, we can just set J=A fixed once, instead of keep update it! Maybe just create a function getFixedJacobian? 161c4762a1bSJed Brown For nonlinear least square, we require x to compute J, keep codes here for future nonlinear least square*/ 162c4762a1bSJed Brown for (m = 0; m < M; ++m) { 163ad540459SPierre Jolivet for (n = 0; n < N; ++n) user->J[m][n] = user->A[m][n]; 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 1669566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, M, user->idm, N, user->idn, (PetscReal *)user->J, INSERT_VALUES)); 1679566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1689566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 169c4762a1bSJed Brown 1709566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); /* not used for linear least square, but keep for future nonlinear least square) */ 171c4762a1bSJed Brown PetscLogFlops(0); /* 0 for linear least square, >0 for nonlinear least square */ 172c4762a1bSJed Brown PetscFunctionReturn(0); 173c4762a1bSJed Brown } 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* ------------------------------------------------------------ */ 176c4762a1bSJed Brown /* Currently fixed matrix, in future may be dynamic for D(x)? */ 177d71ae5a4SJacob Faibussowitsch PetscErrorCode FormDictionaryMatrix(Mat D, AppCtx *user) 178d71ae5a4SJacob Faibussowitsch { 179c4762a1bSJed Brown PetscFunctionBegin; 1809566063dSJacob Faibussowitsch PetscCall(MatSetValues(D, K, user->idk, N, user->idn, (PetscReal *)user->D, INSERT_VALUES)); 1819566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY)); 1829566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY)); 183c4762a1bSJed Brown 184c4762a1bSJed Brown PetscLogFlops(0); /* 0 for fixed dictionary matrix, >0 for varying dictionary matrix */ 185c4762a1bSJed Brown PetscFunctionReturn(0); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* ------------------------------------------------------------ */ 189d71ae5a4SJacob Faibussowitsch PetscErrorCode FormStartingPoint(Vec X) 190d71ae5a4SJacob Faibussowitsch { 191c4762a1bSJed Brown PetscFunctionBegin; 1929566063dSJacob Faibussowitsch PetscCall(VecSet(X, 0.0)); 193c4762a1bSJed Brown PetscFunctionReturn(0); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196c4762a1bSJed Brown /* ---------------------------------------------------------------------- */ 197d71ae5a4SJacob Faibussowitsch PetscErrorCode InitializeUserData(AppCtx *user) 198d71ae5a4SJacob Faibussowitsch { 199*da81f932SPierre Jolivet PetscReal *b = user->b; /* **A=user->A, but we don't know the dimension of A in this way, how to fix? */ 200c4762a1bSJed Brown PetscInt m, n, k; /* loop index for M,N,K dimension. */ 201c4762a1bSJed Brown 202c4762a1bSJed Brown PetscFunctionBegin; 203c4762a1bSJed Brown /* b = A*x while x = [0;0;1;0;0] here*/ 204c4762a1bSJed Brown m = 0; 205c4762a1bSJed Brown b[m++] = 0.28; 206c4762a1bSJed Brown b[m++] = 0.55; 207c4762a1bSJed Brown b[m++] = 0.96; 208c4762a1bSJed Brown 209c4762a1bSJed Brown /* matlab generated random matrix, uniformly distributed in [0,1] with 2 digits accuracy. rng(0); A = rand(M, N); A = round(A*100)/100; 210c4762a1bSJed Brown A = [0.81 0.91 0.28 0.96 0.96 211c4762a1bSJed Brown 0.91 0.63 0.55 0.16 0.49 212c4762a1bSJed Brown 0.13 0.10 0.96 0.97 0.80] 213c4762a1bSJed Brown */ 2149371c9d4SSatish Balay m = 0; 2159371c9d4SSatish Balay n = 0; 2169371c9d4SSatish Balay user->A[m][n++] = 0.81; 2179371c9d4SSatish Balay user->A[m][n++] = 0.91; 2189371c9d4SSatish Balay user->A[m][n++] = 0.28; 2199371c9d4SSatish Balay user->A[m][n++] = 0.96; 2209371c9d4SSatish Balay user->A[m][n++] = 0.96; 2219371c9d4SSatish Balay ++m; 2229371c9d4SSatish Balay n = 0; 2239371c9d4SSatish Balay user->A[m][n++] = 0.91; 2249371c9d4SSatish Balay user->A[m][n++] = 0.63; 2259371c9d4SSatish Balay user->A[m][n++] = 0.55; 2269371c9d4SSatish Balay user->A[m][n++] = 0.16; 2279371c9d4SSatish Balay user->A[m][n++] = 0.49; 2289371c9d4SSatish Balay ++m; 2299371c9d4SSatish Balay n = 0; 2309371c9d4SSatish Balay user->A[m][n++] = 0.13; 2319371c9d4SSatish Balay user->A[m][n++] = 0.10; 2329371c9d4SSatish Balay user->A[m][n++] = 0.96; 2339371c9d4SSatish Balay user->A[m][n++] = 0.97; 2349371c9d4SSatish Balay user->A[m][n++] = 0.80; 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* initialize to 0 */ 237c4762a1bSJed Brown for (k = 0; k < K; k++) { 238ad540459SPierre Jolivet for (n = 0; n < N; n++) user->D[k][n] = 0.0; 239c4762a1bSJed Brown } 240c4762a1bSJed Brown /* Choice I: set D to identity matrix of size N*N for testing */ 241c4762a1bSJed Brown /* for (k=0; k<K; k++) user->D[k][k] = 1.0; */ 242c4762a1bSJed Brown /* Choice II: set D to Backward difference matrix of size (N-1)*N, with zero extended boundary assumption */ 243c4762a1bSJed Brown for (k = 0; k < K; k++) { 244c4762a1bSJed Brown user->D[k][k] = -1.0; 245c4762a1bSJed Brown user->D[k][k + 1] = 1.0; 246c4762a1bSJed Brown } 247c4762a1bSJed Brown 248c4762a1bSJed Brown PetscFunctionReturn(0); 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251c4762a1bSJed Brown /*TEST 252c4762a1bSJed Brown 253c4762a1bSJed Brown build: 254dfd57a17SPierre Jolivet requires: !complex !single !quad !defined(PETSC_USE_64BIT_INDICES) 255c4762a1bSJed Brown 256c4762a1bSJed Brown test: 257c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT 258c4762a1bSJed Brown args: -tao_smonitor -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-6 259c4762a1bSJed Brown 260c4762a1bSJed Brown test: 261c4762a1bSJed Brown suffix: 2 262c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT 2638ebe3e4eSStefano Zampini args: -tao_monitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type l2prox -tao_brgn_regularizer_weight 1e-8 -tao_gatol 1.e-6 -tao_brgn_subsolver_tao_bnk_ksp_converged_reason 264c4762a1bSJed Brown 265c4762a1bSJed Brown test: 266c4762a1bSJed Brown suffix: 3 267c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT 268c4762a1bSJed Brown args: -tao_monitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type l1dict -tao_brgn_regularizer_weight 1e-8 -tao_brgn_l1_smooth_epsilon 1e-6 -tao_gatol 1.e-6 269c4762a1bSJed Brown 270c4762a1bSJed Brown test: 271c4762a1bSJed Brown suffix: 4 272c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT 273c4762a1bSJed Brown args: -tao_monitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type l2pure -tao_brgn_regularizer_weight 1e-8 -tao_gatol 1.e-6 274c4762a1bSJed Brown 275cd1c4666STristan Konolige test: 276cd1c4666STristan Konolige suffix: 5 277cd1c4666STristan Konolige localrunfiles: cs1Data_A_b_xGT 278cd1c4666STristan Konolige args: -tao_monitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type lm -tao_gatol 1.e-6 -tao_brgn_subsolver_tao_type bnls 279cd1c4666STristan Konolige 280c4762a1bSJed Brown TEST*/ 281