xref: /petsc/src/tao/leastsquares/tutorials/tomography.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
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
9c4762a1bSJed Brown      petscsys.h    - sysem 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 /*T
28c4762a1bSJed Brown    Concepts: TAO^Solving a system of nonlinear equations, nonlinear least squares
29c4762a1bSJed Brown    Routines: TaoCreate();
30c4762a1bSJed Brown    Routines: TaoSetType();
31c4762a1bSJed Brown    Routines: TaoSetSeparableObjectiveRoutine();
32c4762a1bSJed Brown    Routines: TaoSetJacobianRoutine();
33c4762a1bSJed Brown    Routines: TaoSetInitialVector();
34c4762a1bSJed Brown    Routines: TaoSetFromOptions();
35c4762a1bSJed Brown    Routines: TaoSetConvergenceHistory(); TaoGetConvergenceHistory();
36c4762a1bSJed Brown    Routines: TaoSolve();
37c4762a1bSJed Brown    Routines: TaoView(); TaoDestroy();
38c4762a1bSJed Brown    Processors: 1
39c4762a1bSJed Brown T*/
40c4762a1bSJed Brown 
41c4762a1bSJed Brown /* User-defined application context */
42c4762a1bSJed Brown typedef struct {
43c4762a1bSJed Brown   /* Working space. linear least square:  res(x) = A*x - b */
44c4762a1bSJed Brown   PetscInt  M,N,K;            /* Problem dimension: A is M*N Matrix, D is K*N Matrix */
45c4762a1bSJed 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 */
46c4762a1bSJed Brown   Vec       b,xGT,xlb,xub;    /* observation b, ground truth xGT, the lower bound and upper bound of x*/
47c4762a1bSJed Brown } AppCtx;
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /* User provided Routines */
50c4762a1bSJed Brown PetscErrorCode InitializeUserData(AppCtx *);
51c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec,AppCtx *);
52c4762a1bSJed Brown PetscErrorCode EvaluateResidual(Tao,Vec,Vec,void *);
53c4762a1bSJed Brown PetscErrorCode EvaluateJacobian(Tao,Vec,Mat,Mat,void *);
54c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerObjectiveAndGradient(Tao,Vec,PetscReal *,Vec,void*);
55c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerHessian(Tao,Vec,Mat,void*);
56c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerHessianProd(Mat,Vec,Vec);
57c4762a1bSJed Brown 
58c4762a1bSJed Brown /*--------------------------------------------------------------------*/
59c4762a1bSJed Brown int main(int argc,char **argv)
60c4762a1bSJed Brown {
61c4762a1bSJed Brown   PetscErrorCode ierr;               /* used to check for functions returning nonzeros */
62c4762a1bSJed Brown   Vec            x,res;              /* solution, function res(x) = A*x-b */
63c4762a1bSJed Brown   Mat            Hreg;               /* regularizer Hessian matrix for user specified regularizer*/
64c4762a1bSJed Brown   Tao            tao;                /* Tao solver context */
65c4762a1bSJed Brown   PetscReal      hist[100],resid[100],v1,v2;
66c4762a1bSJed Brown   PetscInt       lits[100];
67c4762a1bSJed Brown   AppCtx         user;               /* user-defined work context */
68c4762a1bSJed Brown   PetscViewer    fd;   /* used to save result to file */
69c4762a1bSJed Brown   char           resultFile[] = "tomographyResult_x";  /* Debug: change from "tomographyResult_x" to "cs1Result_x" */
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char *)0,help);if (ierr) return ierr;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
74c4762a1bSJed Brown   ierr = TaoCreate(PETSC_COMM_SELF,&tao);CHKERRQ(ierr);
75c4762a1bSJed Brown   ierr = TaoSetType(tao,TAOBRGN);CHKERRQ(ierr);
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* User set application context: A, D matrice, and b vector. */
78c4762a1bSJed Brown   ierr = InitializeUserData(&user);CHKERRQ(ierr);
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Allocate solution vector x,  and function vectors Ax-b, */
81c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,user.N,&x);CHKERRQ(ierr);
82c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,user.M,&res);CHKERRQ(ierr);
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* Set initial guess */
85c4762a1bSJed Brown   ierr = FormStartingPoint(x,&user);CHKERRQ(ierr);
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* Bind x to tao->solution. */
88c4762a1bSJed Brown   ierr = TaoSetInitialVector(tao,x);CHKERRQ(ierr);
89c4762a1bSJed Brown   /* Sets the upper and lower bounds of x */
90c4762a1bSJed Brown   ierr = TaoSetVariableBounds(tao,user.xlb,user.xub);CHKERRQ(ierr);
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   /* Bind user.D to tao->data->D */
93c4762a1bSJed Brown   ierr = TaoBRGNSetDictionaryMatrix(tao,user.D);CHKERRQ(ierr);
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* Set the residual function and Jacobian routines for least squares. */
96c4762a1bSJed Brown   ierr = TaoSetResidualRoutine(tao,res,EvaluateResidual,(void*)&user);CHKERRQ(ierr);
97c4762a1bSJed Brown   /* Jacobian matrix fixed as user.A for Linear least sqaure problem. */
98c4762a1bSJed Brown   ierr = TaoSetJacobianResidualRoutine(tao,user.A,user.A,EvaluateJacobian,(void*)&user);CHKERRQ(ierr);
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* User set the regularizer objective, gradient, and hessian. Set it the same as using l2prox choice, for testing purpose.  */
101c4762a1bSJed Brown   ierr = TaoBRGNSetRegularizerObjectiveAndGradientRoutine(tao,EvaluateRegularizerObjectiveAndGradient,(void*)&user);CHKERRQ(ierr);
102c4762a1bSJed Brown   /* User defined regularizer Hessian setup, here is identiy shell matrix */
103c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_SELF,&Hreg);CHKERRQ(ierr);
104c4762a1bSJed Brown   ierr = MatSetSizes(Hreg,PETSC_DECIDE,PETSC_DECIDE,user.N,user.N);CHKERRQ(ierr);
105c4762a1bSJed Brown   ierr = MatSetType(Hreg,MATSHELL);CHKERRQ(ierr);
106c4762a1bSJed Brown   ierr = MatSetUp(Hreg);CHKERRQ(ierr);
107c4762a1bSJed Brown   ierr = MatShellSetOperation(Hreg,MATOP_MULT,(void (*)(void))EvaluateRegularizerHessianProd);CHKERRQ(ierr);
108c4762a1bSJed Brown   ierr = TaoBRGNSetRegularizerHessianRoutine(tao,Hreg,EvaluateRegularizerHessian,(void*)&user);CHKERRQ(ierr);
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* Check for any TAO command line arguments */
111c4762a1bSJed Brown   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   ierr = TaoSetConvergenceHistory(tao,hist,resid,0,lits,100,PETSC_TRUE);CHKERRQ(ierr);
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* Perform the Solve */
116c4762a1bSJed Brown   ierr = TaoSolve(tao);CHKERRQ(ierr);
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   /* Save x (reconstruction of object) vector to a binary file, which maybe read from Matlab and convert to a 2D image for comparison. */
119c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,resultFile,FILE_MODE_WRITE,&fd);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = VecView(x,fd);CHKERRQ(ierr);
121c4762a1bSJed Brown   ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /* compute the error */
124c4762a1bSJed Brown   ierr = VecAXPY(x,-1,user.xGT);CHKERRQ(ierr);
125c4762a1bSJed Brown   ierr = VecNorm(x,NORM_2,&v1);CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = VecNorm(user.xGT,NORM_2,&v2);CHKERRQ(ierr);
127c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "relative reconstruction error: ||x-xGT||/||xGT|| = %6.4e.\n", (double)(v1/v2));CHKERRQ(ierr);
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* Free TAO data structures */
130c4762a1bSJed Brown   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
131c4762a1bSJed Brown 
132c4762a1bSJed Brown    /* Free PETSc data structures */
133c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
134c4762a1bSJed Brown   ierr = VecDestroy(&res);CHKERRQ(ierr);
135c4762a1bSJed Brown   ierr = MatDestroy(&Hreg);CHKERRQ(ierr);
136c4762a1bSJed Brown   /* Free user data structures */
137c4762a1bSJed Brown   ierr = MatDestroy(&user.A);CHKERRQ(ierr);
138c4762a1bSJed Brown   ierr = MatDestroy(&user.D);CHKERRQ(ierr);
139c4762a1bSJed Brown   ierr = VecDestroy(&user.b);CHKERRQ(ierr);
140c4762a1bSJed Brown   ierr = VecDestroy(&user.xGT);CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr = VecDestroy(&user.xlb);CHKERRQ(ierr);
142c4762a1bSJed Brown   ierr = VecDestroy(&user.xub);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = PetscFinalize();
144c4762a1bSJed Brown   return ierr;
145c4762a1bSJed Brown }
146c4762a1bSJed Brown 
147c4762a1bSJed Brown /*--------------------------------------------------------------------*/
148c4762a1bSJed Brown /* Evaluate residual function A(x)-b in least square problem ||A(x)-b||^2 */
149c4762a1bSJed Brown PetscErrorCode EvaluateResidual(Tao tao,Vec X,Vec F,void *ptr)
150c4762a1bSJed Brown {
151c4762a1bSJed Brown   AppCtx         *user = (AppCtx *)ptr;
152c4762a1bSJed Brown   PetscErrorCode ierr;
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   PetscFunctionBegin;
155c4762a1bSJed Brown   /* Compute Ax - b */
156c4762a1bSJed Brown   ierr = MatMult(user->A,X,F);CHKERRQ(ierr);
157c4762a1bSJed Brown   ierr = VecAXPY(F,-1,user->b);CHKERRQ(ierr);
158ca0c957dSBarry Smith   PetscLogFlops(2.0*user->M*user->N);
159c4762a1bSJed Brown   PetscFunctionReturn(0);
160c4762a1bSJed Brown }
161c4762a1bSJed Brown 
162c4762a1bSJed Brown /*------------------------------------------------------------*/
163c4762a1bSJed Brown PetscErrorCode EvaluateJacobian(Tao tao,Vec X,Mat J,Mat Jpre,void *ptr)
164c4762a1bSJed Brown {
165c4762a1bSJed 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 */
166c4762a1bSJed Brown   PetscFunctionBegin;
167c4762a1bSJed Brown   PetscFunctionReturn(0);
168c4762a1bSJed Brown }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown /* ------------------------------------------------------------ */
171c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerObjectiveAndGradient(Tao tao,Vec X,PetscReal *f_reg,Vec G_reg,void *ptr)
172c4762a1bSJed Brown {
173c4762a1bSJed Brown   PetscErrorCode ierr;
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   PetscFunctionBegin;
176c4762a1bSJed Brown   /* compute regularizer objective = 0.5*x'*x */
177c4762a1bSJed Brown   ierr = VecDot(X,X,f_reg);CHKERRQ(ierr);
178c4762a1bSJed Brown   *f_reg *= 0.5;
179c4762a1bSJed Brown   /* compute regularizer gradient = x */
180c4762a1bSJed Brown   ierr = VecCopy(X,G_reg);CHKERRQ(ierr);
181c4762a1bSJed Brown   PetscFunctionReturn(0);
182c4762a1bSJed Brown }
183c4762a1bSJed Brown 
184c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerHessianProd(Mat Hreg,Vec in,Vec out)
185c4762a1bSJed Brown {
186c4762a1bSJed Brown   PetscErrorCode ierr;
187c4762a1bSJed Brown   PetscFunctionBegin;
188c4762a1bSJed Brown   ierr = VecCopy(in,out);CHKERRQ(ierr);
189c4762a1bSJed Brown   PetscFunctionReturn(0);
190c4762a1bSJed Brown }
191c4762a1bSJed Brown 
192c4762a1bSJed Brown /* ------------------------------------------------------------ */
193c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerHessian(Tao tao,Vec X,Mat Hreg,void *ptr)
194c4762a1bSJed Brown {
195c4762a1bSJed Brown   /* Hessian for regularizer objective = 0.5*x'*x is identity matrix, and is not changing*/
196c4762a1bSJed Brown   PetscFunctionBegin;
197c4762a1bSJed Brown   PetscFunctionReturn(0);
198c4762a1bSJed Brown }
199c4762a1bSJed Brown 
200c4762a1bSJed Brown /* ------------------------------------------------------------ */
201c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec X,AppCtx *user)
202c4762a1bSJed Brown {
203c4762a1bSJed Brown   PetscErrorCode ierr;
204c4762a1bSJed Brown   PetscFunctionBegin;
205c4762a1bSJed Brown   ierr = VecSet(X,0.0);CHKERRQ(ierr);
206c4762a1bSJed Brown   PetscFunctionReturn(0);
207c4762a1bSJed Brown }
208c4762a1bSJed Brown 
209c4762a1bSJed Brown /* ---------------------------------------------------------------------- */
210c4762a1bSJed Brown PetscErrorCode InitializeUserData(AppCtx *user)
211c4762a1bSJed Brown {
212c4762a1bSJed Brown   PetscInt       k,n; /* indices for row and columns of D. */
213c4762a1bSJed 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". */
214c4762a1bSJed Brown   PetscInt       dictChoice = 1; /* choose from 0:identity, 1:gradient1D, 2:gradient2D, 3:DCT etc */
215c4762a1bSJed Brown   PetscViewer    fd;   /* used to load data from file */
216c4762a1bSJed Brown   PetscErrorCode ierr;
217c4762a1bSJed Brown   PetscReal      v;
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   PetscFunctionBegin;
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /*
222c4762a1bSJed Brown   Matrix Vector read and write refer to:
223a17b96a8SKyle Gerard Felker   https://petsc.org/release/src/mat/tutorials/ex10.c
224a17b96a8SKyle Gerard Felker   https://petsc.org/release/src/mat/tutorials/ex12.c
225c4762a1bSJed Brown  */
226c4762a1bSJed Brown   /* Load the A matrix, b vector, and xGT vector from a binary file. */
227c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,dataFile,FILE_MODE_READ,&fd);CHKERRQ(ierr);
228c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&user->A);CHKERRQ(ierr);
229c4762a1bSJed Brown   ierr = MatSetType(user->A,MATSEQAIJ);CHKERRQ(ierr);
230c4762a1bSJed Brown   ierr = MatLoad(user->A,fd);CHKERRQ(ierr);
231c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&user->b);CHKERRQ(ierr);
232c4762a1bSJed Brown   ierr = VecLoad(user->b,fd);CHKERRQ(ierr);
233c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&user->xGT);CHKERRQ(ierr);
234c4762a1bSJed Brown   ierr = VecLoad(user->xGT,fd);CHKERRQ(ierr);
235c4762a1bSJed Brown   ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
236c4762a1bSJed Brown   ierr = VecDuplicate(user->xGT,&(user->xlb));CHKERRQ(ierr);
237c4762a1bSJed Brown   ierr = VecSet(user->xlb,0.0);CHKERRQ(ierr);
238c4762a1bSJed Brown   ierr = VecDuplicate(user->xGT,&(user->xub));CHKERRQ(ierr);
239c4762a1bSJed Brown   ierr = VecSet(user->xub,PETSC_INFINITY);CHKERRQ(ierr);
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   /* Specify the size */
242c4762a1bSJed Brown   ierr = MatGetSize(user->A,&user->M,&user->N);CHKERRQ(ierr);
243c4762a1bSJed Brown 
244c4762a1bSJed 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.
245c4762a1bSJed Brown   if (dictChoice == 0) {
246c4762a1bSJed Brown     user->D = NULL;
247c4762a1bSJed Brown     PetscFunctionReturn(0);
248c4762a1bSJed Brown   }
249c4762a1bSJed Brown   */
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /* Speficy D */
252c4762a1bSJed Brown   /* (1) Specify D Size */
253c4762a1bSJed Brown   switch (dictChoice) {
254c4762a1bSJed Brown     case 0: /* 0:identity */
255c4762a1bSJed Brown       user->K = user->N;
256c4762a1bSJed Brown       break;
257c4762a1bSJed Brown     case 1: /* 1:gradient1D */
258c4762a1bSJed Brown       user->K = user->N-1;
259c4762a1bSJed Brown       break;
260c4762a1bSJed Brown   }
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_SELF,&user->D);CHKERRQ(ierr);
263c4762a1bSJed Brown   ierr = MatSetSizes(user->D,PETSC_DECIDE,PETSC_DECIDE,user->K,user->N);CHKERRQ(ierr);
264c4762a1bSJed Brown   ierr = MatSetFromOptions(user->D);CHKERRQ(ierr);
265c4762a1bSJed Brown   ierr = MatSetUp(user->D);CHKERRQ(ierr);
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /* (2) Specify D Content */
268c4762a1bSJed Brown   switch (dictChoice) {
269c4762a1bSJed Brown     case 0: /* 0:identity */
270c4762a1bSJed Brown       for (k=0; k<user->K; k++) {
271c4762a1bSJed Brown         v = 1.0;
272c4762a1bSJed Brown         ierr = MatSetValues(user->D,1,&k,1,&k,&v,INSERT_VALUES);CHKERRQ(ierr);
273c4762a1bSJed Brown       }
274c4762a1bSJed Brown       break;
275c4762a1bSJed Brown     case 1: /* 1:gradient1D.  [-1, 1, 0,...; 0, -1, 1, 0, ...] */
276c4762a1bSJed Brown       for (k=0; k<user->K; k++) {
277c4762a1bSJed Brown         v = 1.0;
278c4762a1bSJed Brown         n = k+1;
279c4762a1bSJed Brown         ierr = MatSetValues(user->D,1,&k,1,&n,&v,INSERT_VALUES);CHKERRQ(ierr);
280c4762a1bSJed Brown         v = -1.0;
281c4762a1bSJed Brown         ierr = MatSetValues(user->D,1,&k,1,&k,&v,INSERT_VALUES);CHKERRQ(ierr);
282c4762a1bSJed Brown       }
283c4762a1bSJed Brown       break;
284c4762a1bSJed Brown   }
285c4762a1bSJed Brown   ierr = MatAssemblyBegin(user->D,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
286c4762a1bSJed Brown   ierr = MatAssemblyEnd(user->D,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   PetscFunctionReturn(0);
289c4762a1bSJed Brown }
290c4762a1bSJed Brown 
291c4762a1bSJed Brown /*TEST
292c4762a1bSJed Brown 
293c4762a1bSJed Brown    build:
294*dfd57a17SPierre Jolivet       requires: !complex !single !__float128 !defined(PETSC_USE_64BIT_INDICES)
295c4762a1bSJed Brown 
296c4762a1bSJed Brown    test:
297c4762a1bSJed Brown       localrunfiles: tomographyData_A_b_xGT
298c4762a1bSJed 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
299c4762a1bSJed Brown 
300c4762a1bSJed Brown    test:
301c4762a1bSJed Brown       suffix: 2
302c4762a1bSJed Brown       localrunfiles: tomographyData_A_b_xGT
303c4762a1bSJed Brown       args: -tao_monitor -tao_max_it 1000 -tao_brgn_regularization_type l2prox -tao_brgn_regularizer_weight 1e-8 -tao_gatol 1.e-6
304c4762a1bSJed Brown 
305c4762a1bSJed Brown    test:
306c4762a1bSJed Brown       suffix: 3
307c4762a1bSJed Brown       localrunfiles: tomographyData_A_b_xGT
308c4762a1bSJed Brown       args: -tao_monitor -tao_max_it 1000 -tao_brgn_regularization_type user -tao_brgn_regularizer_weight 1e-8 -tao_gatol 1.e-6
309c4762a1bSJed Brown 
310c4762a1bSJed Brown TEST*/
311