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 /*--------------------------------------------------------------------*/ 52*9371c9d4SSatish Balay int main(int argc, char **argv) { 53c4762a1bSJed Brown Vec x, f; /* solution, function f(x) = A*x-b */ 54c4762a1bSJed Brown Mat J, D; /* Jacobian matrix, Transform matrix */ 55c4762a1bSJed Brown Tao tao; /* Tao solver context */ 56c4762a1bSJed Brown PetscInt i; /* iteration information */ 57c4762a1bSJed Brown PetscReal hist[100], resid[100]; 58c4762a1bSJed Brown PetscInt lits[100]; 59c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 60c4762a1bSJed Brown 61327415f7SBarry Smith PetscFunctionBeginUser; 629566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* Allocate solution and vector function vectors */ 659566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &x)); 669566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, M, &f)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Allocate Jacobian and Dictionary matrix. */ 699566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, M, N, NULL, &J)); 709566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, K, N, NULL, &D)); /* XH: TODO: dense -> sparse/dense/shell etc, do it on fly */ 71c4762a1bSJed Brown 72c4762a1bSJed Brown for (i = 0; i < M; i++) user.idm[i] = i; 73c4762a1bSJed Brown for (i = 0; i < N; i++) user.idn[i] = i; 74c4762a1bSJed Brown for (i = 0; i < K; i++) user.idk[i] = i; 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 779566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF, &tao)); 789566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBRGN)); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* User set application context: A, D matrice, and b vector. */ 819566063dSJacob Faibussowitsch PetscCall(InitializeUserData(&user)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* Set initial guess */ 849566063dSJacob Faibussowitsch PetscCall(FormStartingPoint(x)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Fill the content of matrix D from user application Context */ 879566063dSJacob Faibussowitsch PetscCall(FormDictionaryMatrix(D, &user)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* Bind x to tao->solution. */ 909566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 91c4762a1bSJed Brown /* Bind D to tao->data->D */ 929566063dSJacob Faibussowitsch PetscCall(TaoBRGNSetDictionaryMatrix(tao, D)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* Set the function and Jacobian routines. */ 959566063dSJacob Faibussowitsch PetscCall(TaoSetResidualRoutine(tao, f, EvaluateFunction, (void *)&user)); 969566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianResidualRoutine(tao, J, J, EvaluateJacobian, (void *)&user)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* Check for any TAO command line arguments */ 999566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 100c4762a1bSJed Brown 1019566063dSJacob Faibussowitsch PetscCall(TaoSetConvergenceHistory(tao, hist, resid, 0, lits, 100, PETSC_TRUE)); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Perform the Solve */ 1049566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* XH: Debug: View the result, function and Jacobian. */ 1079566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "-------- result x, residual f=A*x-b, and Jacobian=A. -------- \n")); 1089566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_SELF)); 1099566063dSJacob Faibussowitsch PetscCall(VecView(f, PETSC_VIEWER_STDOUT_SELF)); 1109566063dSJacob Faibussowitsch PetscCall(MatView(J, PETSC_VIEWER_STDOUT_SELF)); 1119566063dSJacob Faibussowitsch PetscCall(MatView(D, PETSC_VIEWER_STDOUT_SELF)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Free TAO data structures */ 1149566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* Free PETSc data structures */ 1179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&f)); 1199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 121c4762a1bSJed Brown 1229566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 123b122ec5aSJacob Faibussowitsch return 0; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown /*--------------------------------------------------------------------*/ 127*9371c9d4SSatish Balay PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr) { 128c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 129c4762a1bSJed Brown PetscInt m, n; 130c4762a1bSJed Brown const PetscReal *x; 131c4762a1bSJed Brown PetscReal *b = user->b, *f; 132c4762a1bSJed Brown 133c4762a1bSJed Brown PetscFunctionBegin; 1349566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 1359566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 136c4762a1bSJed Brown 137a5b23f4aSJose 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 */ 138c4762a1bSJed Brown for (m = 0; m < M; m++) { 139c4762a1bSJed Brown f[m] = -b[m]; 140*9371c9d4SSatish Balay for (n = 0; n < N; n++) { f[m] += user->A[m][n] * x[n]; } 141c4762a1bSJed Brown } 1429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 144ca0c957dSBarry Smith PetscLogFlops(2.0 * M * N); 145c4762a1bSJed Brown PetscFunctionReturn(0); 146c4762a1bSJed Brown } 147c4762a1bSJed Brown 148c4762a1bSJed Brown /*------------------------------------------------------------*/ 149c4762a1bSJed Brown /* J[m][n] = df[m]/dx[n] */ 150*9371c9d4SSatish Balay PetscErrorCode EvaluateJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr) { 151c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 152c4762a1bSJed Brown PetscInt m, n; 153c4762a1bSJed Brown const PetscReal *x; 154c4762a1bSJed Brown 155c4762a1bSJed Brown PetscFunctionBegin; 1569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); /* not used for linear least square, but keep for future nonlinear least square) */ 157c4762a1bSJed 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? 158c4762a1bSJed Brown For nonlinear least square, we require x to compute J, keep codes here for future nonlinear least square*/ 159c4762a1bSJed Brown for (m = 0; m < M; ++m) { 160*9371c9d4SSatish Balay for (n = 0; n < N; ++n) { user->J[m][n] = user->A[m][n]; } 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 1639566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, M, user->idm, N, user->idn, (PetscReal *)user->J, INSERT_VALUES)); 1649566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1659566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 166c4762a1bSJed Brown 1679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); /* not used for linear least square, but keep for future nonlinear least square) */ 168c4762a1bSJed Brown PetscLogFlops(0); /* 0 for linear least square, >0 for nonlinear least square */ 169c4762a1bSJed Brown PetscFunctionReturn(0); 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* ------------------------------------------------------------ */ 173c4762a1bSJed Brown /* Currently fixed matrix, in future may be dynamic for D(x)? */ 174*9371c9d4SSatish Balay PetscErrorCode FormDictionaryMatrix(Mat D, AppCtx *user) { 175c4762a1bSJed Brown PetscFunctionBegin; 1769566063dSJacob Faibussowitsch PetscCall(MatSetValues(D, K, user->idk, N, user->idn, (PetscReal *)user->D, INSERT_VALUES)); 1779566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY)); 1789566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown PetscLogFlops(0); /* 0 for fixed dictionary matrix, >0 for varying dictionary matrix */ 181c4762a1bSJed Brown PetscFunctionReturn(0); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* ------------------------------------------------------------ */ 185*9371c9d4SSatish Balay PetscErrorCode FormStartingPoint(Vec X) { 186c4762a1bSJed Brown PetscFunctionBegin; 1879566063dSJacob Faibussowitsch PetscCall(VecSet(X, 0.0)); 188c4762a1bSJed Brown PetscFunctionReturn(0); 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* ---------------------------------------------------------------------- */ 192*9371c9d4SSatish Balay PetscErrorCode InitializeUserData(AppCtx *user) { 193c4762a1bSJed Brown PetscReal *b = user->b; /* **A=user->A, but we don't kown the dimension of A in this way, how to fix? */ 194c4762a1bSJed Brown PetscInt m, n, k; /* loop index for M,N,K dimension. */ 195c4762a1bSJed Brown 196c4762a1bSJed Brown PetscFunctionBegin; 197c4762a1bSJed Brown /* b = A*x while x = [0;0;1;0;0] here*/ 198c4762a1bSJed Brown m = 0; 199c4762a1bSJed Brown b[m++] = 0.28; 200c4762a1bSJed Brown b[m++] = 0.55; 201c4762a1bSJed Brown b[m++] = 0.96; 202c4762a1bSJed Brown 203c4762a1bSJed 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; 204c4762a1bSJed Brown A = [0.81 0.91 0.28 0.96 0.96 205c4762a1bSJed Brown 0.91 0.63 0.55 0.16 0.49 206c4762a1bSJed Brown 0.13 0.10 0.96 0.97 0.80] 207c4762a1bSJed Brown */ 208*9371c9d4SSatish Balay m = 0; 209*9371c9d4SSatish Balay n = 0; 210*9371c9d4SSatish Balay user->A[m][n++] = 0.81; 211*9371c9d4SSatish Balay user->A[m][n++] = 0.91; 212*9371c9d4SSatish Balay user->A[m][n++] = 0.28; 213*9371c9d4SSatish Balay user->A[m][n++] = 0.96; 214*9371c9d4SSatish Balay user->A[m][n++] = 0.96; 215*9371c9d4SSatish Balay ++m; 216*9371c9d4SSatish Balay n = 0; 217*9371c9d4SSatish Balay user->A[m][n++] = 0.91; 218*9371c9d4SSatish Balay user->A[m][n++] = 0.63; 219*9371c9d4SSatish Balay user->A[m][n++] = 0.55; 220*9371c9d4SSatish Balay user->A[m][n++] = 0.16; 221*9371c9d4SSatish Balay user->A[m][n++] = 0.49; 222*9371c9d4SSatish Balay ++m; 223*9371c9d4SSatish Balay n = 0; 224*9371c9d4SSatish Balay user->A[m][n++] = 0.13; 225*9371c9d4SSatish Balay user->A[m][n++] = 0.10; 226*9371c9d4SSatish Balay user->A[m][n++] = 0.96; 227*9371c9d4SSatish Balay user->A[m][n++] = 0.97; 228*9371c9d4SSatish Balay user->A[m][n++] = 0.80; 229c4762a1bSJed Brown 230c4762a1bSJed Brown /* initialize to 0 */ 231c4762a1bSJed Brown for (k = 0; k < K; k++) { 232*9371c9d4SSatish Balay for (n = 0; n < N; n++) { user->D[k][n] = 0.0; } 233c4762a1bSJed Brown } 234c4762a1bSJed Brown /* Choice I: set D to identity matrix of size N*N for testing */ 235c4762a1bSJed Brown /* for (k=0; k<K; k++) user->D[k][k] = 1.0; */ 236c4762a1bSJed Brown /* Choice II: set D to Backward difference matrix of size (N-1)*N, with zero extended boundary assumption */ 237c4762a1bSJed Brown for (k = 0; k < K; k++) { 238c4762a1bSJed Brown user->D[k][k] = -1.0; 239c4762a1bSJed Brown user->D[k][k + 1] = 1.0; 240c4762a1bSJed Brown } 241c4762a1bSJed Brown 242c4762a1bSJed Brown PetscFunctionReturn(0); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown 245c4762a1bSJed Brown /*TEST 246c4762a1bSJed Brown 247c4762a1bSJed Brown build: 248dfd57a17SPierre Jolivet requires: !complex !single !quad !defined(PETSC_USE_64BIT_INDICES) 249c4762a1bSJed Brown 250c4762a1bSJed Brown test: 251c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT 252c4762a1bSJed Brown args: -tao_smonitor -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-6 253c4762a1bSJed Brown 254c4762a1bSJed Brown test: 255c4762a1bSJed Brown suffix: 2 256c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT 2578ebe3e4eSStefano 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 258c4762a1bSJed Brown 259c4762a1bSJed Brown test: 260c4762a1bSJed Brown suffix: 3 261c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT 262c4762a1bSJed 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 263c4762a1bSJed Brown 264c4762a1bSJed Brown test: 265c4762a1bSJed Brown suffix: 4 266c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT 267c4762a1bSJed 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 268c4762a1bSJed Brown 269cd1c4666STristan Konolige test: 270cd1c4666STristan Konolige suffix: 5 271cd1c4666STristan Konolige localrunfiles: cs1Data_A_b_xGT 272cd1c4666STristan 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 273cd1c4666STristan Konolige 274c4762a1bSJed Brown TEST*/ 275