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