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 /*T 27c4762a1bSJed Brown Concepts: TAO^Solving a system of nonlinear equations, nonlinear least squares 28c4762a1bSJed Brown Routines: TaoCreate(); 29c4762a1bSJed Brown Routines: TaoSetType(); 30c4762a1bSJed Brown Routines: TaoSetSeparableObjectiveRoutine(); 31c4762a1bSJed Brown Routines: TaoSetJacobianRoutine(); 32a82e8c82SStefano Zampini Routines: TaoSetSolution(); 33c4762a1bSJed Brown Routines: TaoSetFromOptions(); 34c4762a1bSJed Brown Routines: TaoSetConvergenceHistory(); TaoGetConvergenceHistory(); 35c4762a1bSJed Brown Routines: TaoSolve(); 36c4762a1bSJed Brown Routines: TaoView(); TaoDestroy(); 37c4762a1bSJed Brown Processors: 1 38c4762a1bSJed Brown T*/ 39c4762a1bSJed Brown 40c4762a1bSJed Brown #define M 3 41c4762a1bSJed Brown #define N 5 42c4762a1bSJed Brown #define K 4 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* User-defined application context */ 45c4762a1bSJed Brown typedef struct { 46c4762a1bSJed Brown /* Working space. linear least square: f(x) = A*x - b */ 47c4762a1bSJed Brown PetscReal A[M][N]; /* array of coefficients */ 48c4762a1bSJed Brown PetscReal b[M]; /* array of observations */ 49c4762a1bSJed Brown PetscReal xGT[M]; /* array of ground truth object, which can be used to compare the reconstruction result */ 50c4762a1bSJed Brown PetscReal D[K][N]; /* array of coefficients for 0.5*||Ax-b||^2 + lambda*||D*x||_1 */ 51c4762a1bSJed Brown PetscReal J[M][N]; /* dense jacobian matrix array. For linear least square, J = A. For nonlinear least square, it is different from A */ 52c4762a1bSJed Brown PetscInt idm[M]; /* Matrix row, column indices for jacobian and dictionary */ 53c4762a1bSJed Brown PetscInt idn[N]; 54c4762a1bSJed Brown PetscInt idk[K]; 55c4762a1bSJed Brown } AppCtx; 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* User provided Routines */ 58c4762a1bSJed Brown PetscErrorCode InitializeUserData(AppCtx *); 59c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec); 60c4762a1bSJed Brown PetscErrorCode FormDictionaryMatrix(Mat,AppCtx *); 61c4762a1bSJed Brown PetscErrorCode EvaluateFunction(Tao,Vec,Vec,void *); 62c4762a1bSJed Brown PetscErrorCode EvaluateJacobian(Tao,Vec,Mat,Mat,void *); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /*--------------------------------------------------------------------*/ 65c4762a1bSJed Brown int main(int argc,char **argv) 66c4762a1bSJed Brown { 67c4762a1bSJed Brown PetscErrorCode ierr; /* used to check for functions returning nonzeros */ 68c4762a1bSJed Brown Vec x,f; /* solution, function f(x) = A*x-b */ 69c4762a1bSJed Brown Mat J,D; /* Jacobian matrix, Transform matrix */ 70c4762a1bSJed Brown Tao tao; /* Tao solver context */ 71c4762a1bSJed Brown PetscInt i; /* iteration information */ 72c4762a1bSJed Brown PetscReal hist[100],resid[100]; 73c4762a1bSJed Brown PetscInt lits[100]; 74c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 75c4762a1bSJed Brown 76c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char *)0,help);if (ierr) return ierr; 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Allocate solution and vector function vectors */ 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N,&x)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,M,&f)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* Allocate Jacobian and Dictionary matrix. */ 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,M,N,NULL,&J)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,K,N,NULL,&D)); /* XH: TODO: dense -> sparse/dense/shell etc, do it on fly */ 85c4762a1bSJed Brown 86c4762a1bSJed Brown for (i=0;i<M;i++) user.idm[i] = i; 87c4762a1bSJed Brown for (i=0;i<N;i++) user.idn[i] = i; 88c4762a1bSJed Brown for (i=0;i<K;i++) user.idk[i] = i; 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoCreate(PETSC_COMM_SELF,&tao)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetType(tao,TAOBRGN)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* User set application context: A, D matrice, and b vector. */ 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(InitializeUserData(&user)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Set initial guess */ 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(FormStartingPoint(x)); 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* Fill the content of matrix D from user application Context */ 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(FormDictionaryMatrix(D,&user)); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Bind x to tao->solution. */ 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetSolution(tao,x)); 105c4762a1bSJed Brown /* Bind D to tao->data->D */ 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoBRGNSetDictionaryMatrix(tao,D)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* Set the function and Jacobian routines. */ 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetResidualRoutine(tao,f,EvaluateFunction,(void*)&user)); 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetJacobianResidualRoutine(tao,J,J,EvaluateJacobian,(void*)&user)); 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* Check for any TAO command line arguments */ 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetFromOptions(tao)); 114c4762a1bSJed Brown 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetConvergenceHistory(tao,hist,resid,0,lits,100,PETSC_TRUE)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* Perform the Solve */ 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSolve(tao)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* XH: Debug: View the result, function and Jacobian. */ 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "-------- result x, residual f=A*x-b, and Jacobian=A. -------- \n")); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_SELF)); 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(f,PETSC_VIEWER_STDOUT_SELF)); 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(J,PETSC_VIEWER_STDOUT_SELF)); 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(D,PETSC_VIEWER_STDOUT_SELF)); 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* Free TAO data structures */ 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(TaoDestroy(&tao)); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* Free PETSc data structures */ 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&f)); 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown ierr = PetscFinalize(); 137c4762a1bSJed Brown return ierr; 138c4762a1bSJed Brown } 139c4762a1bSJed Brown 140c4762a1bSJed Brown /*--------------------------------------------------------------------*/ 141c4762a1bSJed Brown PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr) 142c4762a1bSJed Brown { 143c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 144c4762a1bSJed Brown PetscInt m,n; 145c4762a1bSJed Brown const PetscReal *x; 146c4762a1bSJed Brown PetscReal *b=user->b,*f; 147c4762a1bSJed Brown 148c4762a1bSJed Brown PetscFunctionBegin; 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&x)); 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(F,&f)); 151c4762a1bSJed Brown 152a5b23f4aSJose 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 */ 153c4762a1bSJed Brown for (m=0;m<M;m++) { 154c4762a1bSJed Brown f[m] = -b[m]; 155c4762a1bSJed Brown for (n=0;n<N;n++) { 156c4762a1bSJed Brown f[m] += user->A[m][n]*x[n]; 157c4762a1bSJed Brown } 158c4762a1bSJed Brown } 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x)); 160*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(F,&f)); 161ca0c957dSBarry Smith PetscLogFlops(2.0*M*N); 162c4762a1bSJed Brown PetscFunctionReturn(0); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 165c4762a1bSJed Brown /*------------------------------------------------------------*/ 166c4762a1bSJed Brown /* J[m][n] = df[m]/dx[n] */ 167c4762a1bSJed Brown PetscErrorCode EvaluateJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr) 168c4762a1bSJed Brown { 169c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 170c4762a1bSJed Brown PetscInt m,n; 171c4762a1bSJed Brown const PetscReal *x; 172c4762a1bSJed Brown 173c4762a1bSJed Brown PetscFunctionBegin; 174*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&x)); /* not used for linear least square, but keep for future nonlinear least square) */ 175c4762a1bSJed 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? 176c4762a1bSJed Brown For nonlinear least square, we require x to compute J, keep codes here for future nonlinear least square*/ 177c4762a1bSJed Brown for (m=0; m<M; ++m) { 178c4762a1bSJed Brown for (n=0; n<N; ++n) { 179c4762a1bSJed Brown user->J[m][n] = user->A[m][n]; 180c4762a1bSJed Brown } 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(J,M,user->idm,N,user->idn,(PetscReal *)user->J,INSERT_VALUES)); 184*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 185*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 186c4762a1bSJed Brown 187*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x));/* not used for linear least square, but keep for future nonlinear least square) */ 188c4762a1bSJed Brown PetscLogFlops(0); /* 0 for linear least square, >0 for nonlinear least square */ 189c4762a1bSJed Brown PetscFunctionReturn(0); 190c4762a1bSJed Brown } 191c4762a1bSJed Brown 192c4762a1bSJed Brown /* ------------------------------------------------------------ */ 193c4762a1bSJed Brown /* Currently fixed matrix, in future may be dynamic for D(x)? */ 194c4762a1bSJed Brown PetscErrorCode FormDictionaryMatrix(Mat D,AppCtx *user) 195c4762a1bSJed Brown { 196c4762a1bSJed Brown PetscFunctionBegin; 197*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(D,K,user->idk,N,user->idn,(PetscReal *)user->D,INSERT_VALUES)); 198*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(D,MAT_FINAL_ASSEMBLY)); 199*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(D,MAT_FINAL_ASSEMBLY)); 200c4762a1bSJed Brown 201c4762a1bSJed Brown PetscLogFlops(0); /* 0 for fixed dictionary matrix, >0 for varying dictionary matrix */ 202c4762a1bSJed Brown PetscFunctionReturn(0); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* ------------------------------------------------------------ */ 206c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec X) 207c4762a1bSJed Brown { 208c4762a1bSJed Brown PetscFunctionBegin; 209*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(X,0.0)); 210c4762a1bSJed Brown PetscFunctionReturn(0); 211c4762a1bSJed Brown } 212c4762a1bSJed Brown 213c4762a1bSJed Brown /* ---------------------------------------------------------------------- */ 214c4762a1bSJed Brown PetscErrorCode InitializeUserData(AppCtx *user) 215c4762a1bSJed Brown { 216c4762a1bSJed Brown PetscReal *b=user->b; /* **A=user->A, but we don't kown the dimension of A in this way, how to fix? */ 217c4762a1bSJed Brown PetscInt m,n,k; /* loop index for M,N,K dimension. */ 218c4762a1bSJed Brown 219c4762a1bSJed Brown PetscFunctionBegin; 220c4762a1bSJed Brown /* b = A*x while x = [0;0;1;0;0] here*/ 221c4762a1bSJed Brown m = 0; 222c4762a1bSJed Brown b[m++] = 0.28; 223c4762a1bSJed Brown b[m++] = 0.55; 224c4762a1bSJed Brown b[m++] = 0.96; 225c4762a1bSJed Brown 226c4762a1bSJed 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; 227c4762a1bSJed Brown A = [0.81 0.91 0.28 0.96 0.96 228c4762a1bSJed Brown 0.91 0.63 0.55 0.16 0.49 229c4762a1bSJed Brown 0.13 0.10 0.96 0.97 0.80] 230c4762a1bSJed Brown */ 231c4762a1bSJed Brown m=0; n=0; user->A[m][n++] = 0.81; user->A[m][n++] = 0.91; user->A[m][n++] = 0.28; user->A[m][n++] = 0.96; user->A[m][n++] = 0.96; 232c4762a1bSJed Brown ++m; n=0; user->A[m][n++] = 0.91; user->A[m][n++] = 0.63; user->A[m][n++] = 0.55; user->A[m][n++] = 0.16; user->A[m][n++] = 0.49; 233c4762a1bSJed Brown ++m; n=0; user->A[m][n++] = 0.13; user->A[m][n++] = 0.10; user->A[m][n++] = 0.96; user->A[m][n++] = 0.97; user->A[m][n++] = 0.80; 234c4762a1bSJed Brown 235c4762a1bSJed Brown /* initialize to 0 */ 236c4762a1bSJed Brown for (k=0; k<K; k++) { 237c4762a1bSJed Brown for (n=0; n<N; n++) { 238c4762a1bSJed Brown user->D[k][n] = 0.0; 239c4762a1bSJed Brown } 240c4762a1bSJed Brown } 241c4762a1bSJed Brown /* Choice I: set D to identity matrix of size N*N for testing */ 242c4762a1bSJed Brown /* for (k=0; k<K; k++) user->D[k][k] = 1.0; */ 243c4762a1bSJed Brown /* Choice II: set D to Backward difference matrix of size (N-1)*N, with zero extended boundary assumption */ 244c4762a1bSJed Brown for (k=0;k<K;k++) { 245c4762a1bSJed Brown user->D[k][k] = -1.0; 246c4762a1bSJed Brown user->D[k][k+1] = 1.0; 247c4762a1bSJed Brown } 248c4762a1bSJed Brown 249c4762a1bSJed Brown PetscFunctionReturn(0); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252c4762a1bSJed Brown /*TEST 253c4762a1bSJed Brown 254c4762a1bSJed Brown build: 255dfd57a17SPierre Jolivet requires: !complex !single !quad !defined(PETSC_USE_64BIT_INDICES) 256c4762a1bSJed Brown 257c4762a1bSJed Brown test: 258c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT 259c4762a1bSJed Brown args: -tao_smonitor -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-6 260c4762a1bSJed Brown 261c4762a1bSJed Brown test: 262c4762a1bSJed Brown suffix: 2 263c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT 2648ebe3e4eSStefano 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 265c4762a1bSJed Brown 266c4762a1bSJed Brown test: 267c4762a1bSJed Brown suffix: 3 268c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT 269c4762a1bSJed 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 270c4762a1bSJed Brown 271c4762a1bSJed Brown test: 272c4762a1bSJed Brown suffix: 4 273c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT 274c4762a1bSJed 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 275c4762a1bSJed Brown 276cd1c4666STristan Konolige test: 277cd1c4666STristan Konolige suffix: 5 278cd1c4666STristan Konolige localrunfiles: cs1Data_A_b_xGT 279cd1c4666STristan 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 280cd1c4666STristan Konolige 281c4762a1bSJed Brown TEST*/ 282