1c4762a1bSJed Brown /* XH: 2c4762a1bSJed Brown Todo: add cs1f.F90 and adjust makefile. 3c4762a1bSJed Brown Todo: maybe provide code template to generate 1D/2D/3D gradient, DCT transform matrix for D etc. 4c4762a1bSJed Brown */ 5c4762a1bSJed Brown /* 6c4762a1bSJed Brown Include "petsctao.h" so that we can use TAO solvers. Note that this 7c4762a1bSJed Brown file automatically includes libraries such as: 8c4762a1bSJed Brown petsc.h - base PETSc routines petscvec.h - vectors 9a5b23f4aSJose E. Roman petscsys.h - system routines petscmat.h - matrices 10c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 11c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 12c4762a1bSJed Brown 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown 15c4762a1bSJed Brown #include <petsctao.h> 16c4762a1bSJed Brown 17c4762a1bSJed Brown /* 18c4762a1bSJed Brown Description: BRGN tomography reconstruction example . 19c4762a1bSJed Brown 0.5*||Ax-b||^2 + lambda*g(x) 20c4762a1bSJed Brown Reference: None 21c4762a1bSJed Brown */ 22c4762a1bSJed Brown 23c4762a1bSJed Brown static char help[] = "Finds the least-squares solution to the under constraint linear model Ax = b, with regularizer. \n\ 24c4762a1bSJed Brown A is a M*N real matrix (M<N), x is sparse. A good regularizer is an L1 regularizer. \n\ 25c4762a1bSJed 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\ 26c4762a1bSJed 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"; 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* User-defined application context */ 29c4762a1bSJed Brown typedef struct { 30c4762a1bSJed Brown /* Working space. linear least square: res(x) = A*x - b */ 31c4762a1bSJed Brown PetscInt M,N,K; /* Problem dimension: A is M*N Matrix, D is K*N Matrix */ 32c4762a1bSJed Brown Mat A,D; /* Coefficients, Dictionary Transform of size M*N and K*N respectively. For linear least square, Jacobian Matrix J = A. For nonlinear least square, it is different from A */ 33c4762a1bSJed Brown Vec b,xGT,xlb,xub; /* observation b, ground truth xGT, the lower bound and upper bound of x*/ 34c4762a1bSJed Brown } AppCtx; 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* User provided Routines */ 37c4762a1bSJed Brown PetscErrorCode InitializeUserData(AppCtx *); 38c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec,AppCtx *); 39c4762a1bSJed Brown PetscErrorCode EvaluateResidual(Tao,Vec,Vec,void *); 40c4762a1bSJed Brown PetscErrorCode EvaluateJacobian(Tao,Vec,Mat,Mat,void *); 41c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerObjectiveAndGradient(Tao,Vec,PetscReal *,Vec,void*); 42c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerHessian(Tao,Vec,Mat,void*); 43c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerHessianProd(Mat,Vec,Vec); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /*--------------------------------------------------------------------*/ 46c4762a1bSJed Brown int main(int argc,char **argv) 47c4762a1bSJed Brown { 48c4762a1bSJed Brown Vec x,res; /* solution, function res(x) = A*x-b */ 49c4762a1bSJed Brown Mat Hreg; /* regularizer Hessian matrix for user specified regularizer*/ 50c4762a1bSJed Brown Tao tao; /* Tao solver context */ 51c4762a1bSJed Brown PetscReal hist[100],resid[100],v1,v2; 52c4762a1bSJed Brown PetscInt lits[100]; 53c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 54c4762a1bSJed Brown PetscViewer fd; /* used to save result to file */ 55c4762a1bSJed Brown char resultFile[] = "tomographyResult_x"; /* Debug: change from "tomographyResult_x" to "cs1Result_x" */ 56c4762a1bSJed Brown 57*327415f7SBarry Smith PetscFunctionBeginUser; 589566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char *)0,help)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 619566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF,&tao)); 629566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAOBRGN)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* User set application context: A, D matrice, and b vector. */ 659566063dSJacob Faibussowitsch PetscCall(InitializeUserData(&user)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* Allocate solution vector x, and function vectors Ax-b, */ 689566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,user.N,&x)); 699566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,user.M,&res)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* Set initial guess */ 729566063dSJacob Faibussowitsch PetscCall(FormStartingPoint(x,&user)); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* Bind x to tao->solution. */ 759566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,x)); 76c4762a1bSJed Brown /* Sets the upper and lower bounds of x */ 779566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(tao,user.xlb,user.xub)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Bind user.D to tao->data->D */ 809566063dSJacob Faibussowitsch PetscCall(TaoBRGNSetDictionaryMatrix(tao,user.D)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* Set the residual function and Jacobian routines for least squares. */ 839566063dSJacob Faibussowitsch PetscCall(TaoSetResidualRoutine(tao,res,EvaluateResidual,(void*)&user)); 84a5b23f4aSJose E. Roman /* Jacobian matrix fixed as user.A for Linear least square problem. */ 859566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianResidualRoutine(tao,user.A,user.A,EvaluateJacobian,(void*)&user)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* User set the regularizer objective, gradient, and hessian. Set it the same as using l2prox choice, for testing purpose. */ 889566063dSJacob Faibussowitsch PetscCall(TaoBRGNSetRegularizerObjectiveAndGradientRoutine(tao,EvaluateRegularizerObjectiveAndGradient,(void*)&user)); 89c4762a1bSJed Brown /* User defined regularizer Hessian setup, here is identiy shell matrix */ 909566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&Hreg)); 919566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Hreg,PETSC_DECIDE,PETSC_DECIDE,user.N,user.N)); 929566063dSJacob Faibussowitsch PetscCall(MatSetType(Hreg,MATSHELL)); 939566063dSJacob Faibussowitsch PetscCall(MatSetUp(Hreg)); 949566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(Hreg,MATOP_MULT,(void (*)(void))EvaluateRegularizerHessianProd)); 959566063dSJacob Faibussowitsch PetscCall(TaoBRGNSetRegularizerHessianRoutine(tao,Hreg,EvaluateRegularizerHessian,(void*)&user)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Check for any TAO command line arguments */ 989566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 99c4762a1bSJed Brown 1009566063dSJacob Faibussowitsch PetscCall(TaoSetConvergenceHistory(tao,hist,resid,0,lits,100,PETSC_TRUE)); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* Perform the Solve */ 1039566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Save x (reconstruction of object) vector to a binary file, which maybe read from Matlab and convert to a 2D image for comparison. */ 1069566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF,resultFile,FILE_MODE_WRITE,&fd)); 1079566063dSJacob Faibussowitsch PetscCall(VecView(x,fd)); 1089566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* compute the error */ 1119566063dSJacob Faibussowitsch PetscCall(VecAXPY(x,-1,user.xGT)); 1129566063dSJacob Faibussowitsch PetscCall(VecNorm(x,NORM_2,&v1)); 1139566063dSJacob Faibussowitsch PetscCall(VecNorm(user.xGT,NORM_2,&v2)); 1149566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "relative reconstruction error: ||x-xGT||/||xGT|| = %6.4e.\n", (double)(v1/v2))); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* Free TAO data structures */ 1179566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* Free PETSc data structures */ 1209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&res)); 1229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Hreg)); 123c4762a1bSJed Brown /* Free user data structures */ 1249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.A)); 1259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.D)); 1269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.b)); 1279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.xGT)); 1289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.xlb)); 1299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.xub)); 1309566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 131b122ec5aSJacob Faibussowitsch return 0; 132c4762a1bSJed Brown } 133c4762a1bSJed Brown 134c4762a1bSJed Brown /*--------------------------------------------------------------------*/ 135c4762a1bSJed Brown /* Evaluate residual function A(x)-b in least square problem ||A(x)-b||^2 */ 136c4762a1bSJed Brown PetscErrorCode EvaluateResidual(Tao tao,Vec X,Vec F,void *ptr) 137c4762a1bSJed Brown { 138c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 139c4762a1bSJed Brown 140c4762a1bSJed Brown PetscFunctionBegin; 141c4762a1bSJed Brown /* Compute Ax - b */ 1429566063dSJacob Faibussowitsch PetscCall(MatMult(user->A,X,F)); 1439566063dSJacob Faibussowitsch PetscCall(VecAXPY(F,-1,user->b)); 144ca0c957dSBarry Smith PetscLogFlops(2.0*user->M*user->N); 145c4762a1bSJed Brown PetscFunctionReturn(0); 146c4762a1bSJed Brown } 147c4762a1bSJed Brown 148c4762a1bSJed Brown /*------------------------------------------------------------*/ 149c4762a1bSJed Brown PetscErrorCode EvaluateJacobian(Tao tao,Vec X,Mat J,Mat Jpre,void *ptr) 150c4762a1bSJed Brown { 151c4762a1bSJed Brown /* Jacobian is not changing here, so use a empty dummy function here. J[m][n] = df[m]/dx[n] = A[m][n] for linear least square */ 152c4762a1bSJed Brown PetscFunctionBegin; 153c4762a1bSJed Brown PetscFunctionReturn(0); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* ------------------------------------------------------------ */ 157c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerObjectiveAndGradient(Tao tao,Vec X,PetscReal *f_reg,Vec G_reg,void *ptr) 158c4762a1bSJed Brown { 159c4762a1bSJed Brown PetscFunctionBegin; 160c4762a1bSJed Brown /* compute regularizer objective = 0.5*x'*x */ 1619566063dSJacob Faibussowitsch PetscCall(VecDot(X,X,f_reg)); 162c4762a1bSJed Brown *f_reg *= 0.5; 163c4762a1bSJed Brown /* compute regularizer gradient = x */ 1649566063dSJacob Faibussowitsch PetscCall(VecCopy(X,G_reg)); 165c4762a1bSJed Brown PetscFunctionReturn(0); 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 168c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerHessianProd(Mat Hreg,Vec in,Vec out) 169c4762a1bSJed Brown { 170c4762a1bSJed Brown PetscFunctionBegin; 1719566063dSJacob Faibussowitsch PetscCall(VecCopy(in,out)); 172c4762a1bSJed Brown PetscFunctionReturn(0); 173c4762a1bSJed Brown } 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* ------------------------------------------------------------ */ 176c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerHessian(Tao tao,Vec X,Mat Hreg,void *ptr) 177c4762a1bSJed Brown { 178c4762a1bSJed Brown /* Hessian for regularizer objective = 0.5*x'*x is identity matrix, and is not changing*/ 179c4762a1bSJed Brown PetscFunctionBegin; 180c4762a1bSJed Brown PetscFunctionReturn(0); 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183c4762a1bSJed Brown /* ------------------------------------------------------------ */ 184c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec X,AppCtx *user) 185c4762a1bSJed Brown { 186c4762a1bSJed Brown PetscFunctionBegin; 1879566063dSJacob Faibussowitsch PetscCall(VecSet(X,0.0)); 188c4762a1bSJed Brown PetscFunctionReturn(0); 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* ---------------------------------------------------------------------- */ 192c4762a1bSJed Brown PetscErrorCode InitializeUserData(AppCtx *user) 193c4762a1bSJed Brown { 194c4762a1bSJed Brown PetscInt k,n; /* indices for row and columns of D. */ 195c4762a1bSJed Brown char dataFile[] = "tomographyData_A_b_xGT"; /* Matrix A and vectors b, xGT(ground truth) binary files generated by Matlab. Debug: change from "tomographyData_A_b_xGT" to "cs1Data_A_b_xGT". */ 196c4762a1bSJed Brown PetscInt dictChoice = 1; /* choose from 0:identity, 1:gradient1D, 2:gradient2D, 3:DCT etc */ 197c4762a1bSJed Brown PetscViewer fd; /* used to load data from file */ 198c4762a1bSJed Brown PetscReal v; 199c4762a1bSJed Brown 200c4762a1bSJed Brown PetscFunctionBegin; 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* 203c4762a1bSJed Brown Matrix Vector read and write refer to: 204a17b96a8SKyle Gerard Felker https://petsc.org/release/src/mat/tutorials/ex10.c 205a17b96a8SKyle Gerard Felker https://petsc.org/release/src/mat/tutorials/ex12.c 206c4762a1bSJed Brown */ 207c4762a1bSJed Brown /* Load the A matrix, b vector, and xGT vector from a binary file. */ 2089566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,dataFile,FILE_MODE_READ,&fd)); 2099566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user->A)); 2109566063dSJacob Faibussowitsch PetscCall(MatSetType(user->A,MATSEQAIJ)); 2119566063dSJacob Faibussowitsch PetscCall(MatLoad(user->A,fd)); 2129566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user->b)); 2139566063dSJacob Faibussowitsch PetscCall(VecLoad(user->b,fd)); 2149566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user->xGT)); 2159566063dSJacob Faibussowitsch PetscCall(VecLoad(user->xGT,fd)); 2169566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 2179566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->xGT,&(user->xlb))); 2189566063dSJacob Faibussowitsch PetscCall(VecSet(user->xlb,0.0)); 2199566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->xGT,&(user->xub))); 2209566063dSJacob Faibussowitsch PetscCall(VecSet(user->xub,PETSC_INFINITY)); 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* Specify the size */ 2239566063dSJacob Faibussowitsch PetscCall(MatGetSize(user->A,&user->M,&user->N)); 224c4762a1bSJed Brown 225c4762a1bSJed Brown /* shortcut, when D is identity matrix, we may just specify it as NULL, and brgn will treat D*x as x without actually computing D*x. 226c4762a1bSJed Brown if (dictChoice == 0) { 227c4762a1bSJed Brown user->D = NULL; 228c4762a1bSJed Brown PetscFunctionReturn(0); 229c4762a1bSJed Brown } 230c4762a1bSJed Brown */ 231c4762a1bSJed Brown 232c4762a1bSJed Brown /* Speficy D */ 233c4762a1bSJed Brown /* (1) Specify D Size */ 234c4762a1bSJed Brown switch (dictChoice) { 235c4762a1bSJed Brown case 0: /* 0:identity */ 236c4762a1bSJed Brown user->K = user->N; 237c4762a1bSJed Brown break; 238c4762a1bSJed Brown case 1: /* 1:gradient1D */ 239c4762a1bSJed Brown user->K = user->N-1; 240c4762a1bSJed Brown break; 241c4762a1bSJed Brown } 242c4762a1bSJed Brown 2439566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&user->D)); 2449566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->D,PETSC_DECIDE,PETSC_DECIDE,user->K,user->N)); 2459566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->D)); 2469566063dSJacob Faibussowitsch PetscCall(MatSetUp(user->D)); 247c4762a1bSJed Brown 248c4762a1bSJed Brown /* (2) Specify D Content */ 249c4762a1bSJed Brown switch (dictChoice) { 250c4762a1bSJed Brown case 0: /* 0:identity */ 251c4762a1bSJed Brown for (k=0; k<user->K; k++) { 252c4762a1bSJed Brown v = 1.0; 2539566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->D,1,&k,1,&k,&v,INSERT_VALUES)); 254c4762a1bSJed Brown } 255c4762a1bSJed Brown break; 256c4762a1bSJed Brown case 1: /* 1:gradient1D. [-1, 1, 0,...; 0, -1, 1, 0, ...] */ 257c4762a1bSJed Brown for (k=0; k<user->K; k++) { 258c4762a1bSJed Brown v = 1.0; 259c4762a1bSJed Brown n = k+1; 2609566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->D,1,&k,1,&n,&v,INSERT_VALUES)); 261c4762a1bSJed Brown v = -1.0; 2629566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->D,1,&k,1,&k,&v,INSERT_VALUES)); 263c4762a1bSJed Brown } 264c4762a1bSJed Brown break; 265c4762a1bSJed Brown } 2669566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->D,MAT_FINAL_ASSEMBLY)); 2679566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->D,MAT_FINAL_ASSEMBLY)); 268c4762a1bSJed Brown 269c4762a1bSJed Brown PetscFunctionReturn(0); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 272c4762a1bSJed Brown /*TEST 273c4762a1bSJed Brown 274c4762a1bSJed Brown build: 275dfd57a17SPierre Jolivet requires: !complex !single !__float128 !defined(PETSC_USE_64BIT_INDICES) 276c4762a1bSJed Brown 277c4762a1bSJed Brown test: 278c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 279c4762a1bSJed Brown args: -tao_max_it 1000 -tao_brgn_regularization_type l1dict -tao_brgn_regularizer_weight 1e-8 -tao_brgn_l1_smooth_epsilon 1e-6 -tao_gatol 1.e-8 280c4762a1bSJed Brown 281c4762a1bSJed Brown test: 282c4762a1bSJed Brown suffix: 2 283c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 284c4762a1bSJed Brown args: -tao_monitor -tao_max_it 1000 -tao_brgn_regularization_type l2prox -tao_brgn_regularizer_weight 1e-8 -tao_gatol 1.e-6 285c4762a1bSJed Brown 286c4762a1bSJed Brown test: 287c4762a1bSJed Brown suffix: 3 288c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 289c4762a1bSJed Brown args: -tao_monitor -tao_max_it 1000 -tao_brgn_regularization_type user -tao_brgn_regularizer_weight 1e-8 -tao_gatol 1.e-6 290c4762a1bSJed Brown 291c4762a1bSJed Brown TEST*/ 292