xref: /petsc/src/tao/leastsquares/tutorials/tomography.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 /*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();
33a82e8c82SStefano Zampini    Routines: TaoSetSolution();
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   Vec            x,res;              /* solution, function res(x) = A*x-b */
62c4762a1bSJed Brown   Mat            Hreg;               /* regularizer Hessian matrix for user specified regularizer*/
63c4762a1bSJed Brown   Tao            tao;                /* Tao solver context */
64c4762a1bSJed Brown   PetscReal      hist[100],resid[100],v1,v2;
65c4762a1bSJed Brown   PetscInt       lits[100];
66c4762a1bSJed Brown   AppCtx         user;               /* user-defined work context */
67c4762a1bSJed Brown   PetscViewer    fd;   /* used to save result to file */
68c4762a1bSJed Brown   char           resultFile[] = "tomographyResult_x";  /* Debug: change from "tomographyResult_x" to "cs1Result_x" */
69c4762a1bSJed Brown 
70*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char *)0,help));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
735f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_SELF,&tao));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOBRGN));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* User set application context: A, D matrice, and b vector. */
775f80ce2aSJacob Faibussowitsch   CHKERRQ(InitializeUserData(&user));
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* Allocate solution vector x,  and function vectors Ax-b, */
805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,user.N,&x));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,user.M,&res));
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* Set initial guess */
845f80ce2aSJacob Faibussowitsch   CHKERRQ(FormStartingPoint(x,&user));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* Bind x to tao->solution. */
875f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,x));
88c4762a1bSJed Brown   /* Sets the upper and lower bounds of x */
895f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetVariableBounds(tao,user.xlb,user.xub));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* Bind user.D to tao->data->D */
925f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoBRGNSetDictionaryMatrix(tao,user.D));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /* Set the residual function and Jacobian routines for least squares. */
955f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetResidualRoutine(tao,res,EvaluateResidual,(void*)&user));
96a5b23f4aSJose E. Roman   /* Jacobian matrix fixed as user.A for Linear least square problem. */
975f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetJacobianResidualRoutine(tao,user.A,user.A,EvaluateJacobian,(void*)&user));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* User set the regularizer objective, gradient, and hessian. Set it the same as using l2prox choice, for testing purpose.  */
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoBRGNSetRegularizerObjectiveAndGradientRoutine(tao,EvaluateRegularizerObjectiveAndGradient,(void*)&user));
101c4762a1bSJed Brown   /* User defined regularizer Hessian setup, here is identiy shell matrix */
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&Hreg));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(Hreg,PETSC_DECIDE,PETSC_DECIDE,user.N,user.N));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(Hreg,MATSHELL));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(Hreg));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(Hreg,MATOP_MULT,(void (*)(void))EvaluateRegularizerHessianProd));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoBRGNSetRegularizerHessianRoutine(tao,Hreg,EvaluateRegularizerHessian,(void*)&user));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* Check for any TAO command line arguments */
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
111c4762a1bSJed Brown 
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetConvergenceHistory(tao,hist,resid,0,lits,100,PETSC_TRUE));
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /* Perform the Solve */
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /* Save x (reconstruction of object) vector to a binary file, which maybe read from Matlab and convert to a 2D image for comparison. */
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,resultFile,FILE_MODE_WRITE,&fd));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,fd));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   /* compute the error */
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(x,-1,user.xGT));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(x,NORM_2,&v1));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(user.xGT,NORM_2,&v2));
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "relative reconstruction error: ||x-xGT||/||xGT|| = %6.4e.\n", (double)(v1/v2)));
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   /* Free TAO data structures */
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
130c4762a1bSJed Brown 
131c4762a1bSJed Brown    /* Free PETSc data structures */
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&res));
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Hreg));
135c4762a1bSJed Brown   /* Free user data structures */
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.A));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.D));
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.b));
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.xGT));
1405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.xlb));
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.xub));
142*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
143*b122ec5aSJacob Faibussowitsch   return 0;
144c4762a1bSJed Brown }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown /*--------------------------------------------------------------------*/
147c4762a1bSJed Brown /* Evaluate residual function A(x)-b in least square problem ||A(x)-b||^2 */
148c4762a1bSJed Brown PetscErrorCode EvaluateResidual(Tao tao,Vec X,Vec F,void *ptr)
149c4762a1bSJed Brown {
150c4762a1bSJed Brown   AppCtx         *user = (AppCtx *)ptr;
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   PetscFunctionBegin;
153c4762a1bSJed Brown   /* Compute Ax - b */
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->A,X,F));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(F,-1,user->b));
156ca0c957dSBarry Smith   PetscLogFlops(2.0*user->M*user->N);
157c4762a1bSJed Brown   PetscFunctionReturn(0);
158c4762a1bSJed Brown }
159c4762a1bSJed Brown 
160c4762a1bSJed Brown /*------------------------------------------------------------*/
161c4762a1bSJed Brown PetscErrorCode EvaluateJacobian(Tao tao,Vec X,Mat J,Mat Jpre,void *ptr)
162c4762a1bSJed Brown {
163c4762a1bSJed 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 */
164c4762a1bSJed Brown   PetscFunctionBegin;
165c4762a1bSJed Brown   PetscFunctionReturn(0);
166c4762a1bSJed Brown }
167c4762a1bSJed Brown 
168c4762a1bSJed Brown /* ------------------------------------------------------------ */
169c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerObjectiveAndGradient(Tao tao,Vec X,PetscReal *f_reg,Vec G_reg,void *ptr)
170c4762a1bSJed Brown {
171c4762a1bSJed Brown   PetscFunctionBegin;
172c4762a1bSJed Brown   /* compute regularizer objective = 0.5*x'*x */
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDot(X,X,f_reg));
174c4762a1bSJed Brown   *f_reg *= 0.5;
175c4762a1bSJed Brown   /* compute regularizer gradient = x */
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(X,G_reg));
177c4762a1bSJed Brown   PetscFunctionReturn(0);
178c4762a1bSJed Brown }
179c4762a1bSJed Brown 
180c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerHessianProd(Mat Hreg,Vec in,Vec out)
181c4762a1bSJed Brown {
182c4762a1bSJed Brown   PetscFunctionBegin;
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(in,out));
184c4762a1bSJed Brown   PetscFunctionReturn(0);
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown /* ------------------------------------------------------------ */
188c4762a1bSJed Brown PetscErrorCode EvaluateRegularizerHessian(Tao tao,Vec X,Mat Hreg,void *ptr)
189c4762a1bSJed Brown {
190c4762a1bSJed Brown   /* Hessian for regularizer objective = 0.5*x'*x is identity matrix, and is not changing*/
191c4762a1bSJed Brown   PetscFunctionBegin;
192c4762a1bSJed Brown   PetscFunctionReturn(0);
193c4762a1bSJed Brown }
194c4762a1bSJed Brown 
195c4762a1bSJed Brown /* ------------------------------------------------------------ */
196c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec X,AppCtx *user)
197c4762a1bSJed Brown {
198c4762a1bSJed Brown   PetscFunctionBegin;
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(X,0.0));
200c4762a1bSJed Brown   PetscFunctionReturn(0);
201c4762a1bSJed Brown }
202c4762a1bSJed Brown 
203c4762a1bSJed Brown /* ---------------------------------------------------------------------- */
204c4762a1bSJed Brown PetscErrorCode InitializeUserData(AppCtx *user)
205c4762a1bSJed Brown {
206c4762a1bSJed Brown   PetscInt       k,n; /* indices for row and columns of D. */
207c4762a1bSJed 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". */
208c4762a1bSJed Brown   PetscInt       dictChoice = 1; /* choose from 0:identity, 1:gradient1D, 2:gradient2D, 3:DCT etc */
209c4762a1bSJed Brown   PetscViewer    fd;   /* used to load data from file */
210c4762a1bSJed Brown   PetscReal      v;
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   PetscFunctionBegin;
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   /*
215c4762a1bSJed Brown   Matrix Vector read and write refer to:
216a17b96a8SKyle Gerard Felker   https://petsc.org/release/src/mat/tutorials/ex10.c
217a17b96a8SKyle Gerard Felker   https://petsc.org/release/src/mat/tutorials/ex12.c
218c4762a1bSJed Brown  */
219c4762a1bSJed Brown   /* Load the A matrix, b vector, and xGT vector from a binary file. */
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,dataFile,FILE_MODE_READ,&fd));
2215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user->A));
2225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(user->A,MATSEQAIJ));
2235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(user->A,fd));
2245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&user->b));
2255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecLoad(user->b,fd));
2265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&user->xGT));
2275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecLoad(user->xGT,fd));
2285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
2295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user->xGT,&(user->xlb)));
2305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(user->xlb,0.0));
2315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user->xGT,&(user->xub)));
2325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(user->xub,PETSC_INFINITY));
233c4762a1bSJed Brown 
234c4762a1bSJed Brown   /* Specify the size */
2355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(user->A,&user->M,&user->N));
236c4762a1bSJed Brown 
237c4762a1bSJed 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.
238c4762a1bSJed Brown   if (dictChoice == 0) {
239c4762a1bSJed Brown     user->D = NULL;
240c4762a1bSJed Brown     PetscFunctionReturn(0);
241c4762a1bSJed Brown   }
242c4762a1bSJed Brown   */
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   /* Speficy D */
245c4762a1bSJed Brown   /* (1) Specify D Size */
246c4762a1bSJed Brown   switch (dictChoice) {
247c4762a1bSJed Brown     case 0: /* 0:identity */
248c4762a1bSJed Brown       user->K = user->N;
249c4762a1bSJed Brown       break;
250c4762a1bSJed Brown     case 1: /* 1:gradient1D */
251c4762a1bSJed Brown       user->K = user->N-1;
252c4762a1bSJed Brown       break;
253c4762a1bSJed Brown   }
254c4762a1bSJed Brown 
2555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&user->D));
2565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user->D,PETSC_DECIDE,PETSC_DECIDE,user->K,user->N));
2575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user->D));
2585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(user->D));
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   /* (2) Specify D Content */
261c4762a1bSJed Brown   switch (dictChoice) {
262c4762a1bSJed Brown     case 0: /* 0:identity */
263c4762a1bSJed Brown       for (k=0; k<user->K; k++) {
264c4762a1bSJed Brown         v = 1.0;
2655f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(user->D,1,&k,1,&k,&v,INSERT_VALUES));
266c4762a1bSJed Brown       }
267c4762a1bSJed Brown       break;
268c4762a1bSJed Brown     case 1: /* 1:gradient1D.  [-1, 1, 0,...; 0, -1, 1, 0, ...] */
269c4762a1bSJed Brown       for (k=0; k<user->K; k++) {
270c4762a1bSJed Brown         v = 1.0;
271c4762a1bSJed Brown         n = k+1;
2725f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(user->D,1,&k,1,&n,&v,INSERT_VALUES));
273c4762a1bSJed Brown         v = -1.0;
2745f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(user->D,1,&k,1,&k,&v,INSERT_VALUES));
275c4762a1bSJed Brown       }
276c4762a1bSJed Brown       break;
277c4762a1bSJed Brown   }
2785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(user->D,MAT_FINAL_ASSEMBLY));
2795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(user->D,MAT_FINAL_ASSEMBLY));
280c4762a1bSJed Brown 
281c4762a1bSJed Brown   PetscFunctionReturn(0);
282c4762a1bSJed Brown }
283c4762a1bSJed Brown 
284c4762a1bSJed Brown /*TEST
285c4762a1bSJed Brown 
286c4762a1bSJed Brown    build:
287dfd57a17SPierre Jolivet       requires: !complex !single !__float128 !defined(PETSC_USE_64BIT_INDICES)
288c4762a1bSJed Brown 
289c4762a1bSJed Brown    test:
290c4762a1bSJed Brown       localrunfiles: tomographyData_A_b_xGT
291c4762a1bSJed 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
292c4762a1bSJed Brown 
293c4762a1bSJed Brown    test:
294c4762a1bSJed Brown       suffix: 2
295c4762a1bSJed Brown       localrunfiles: tomographyData_A_b_xGT
296c4762a1bSJed Brown       args: -tao_monitor -tao_max_it 1000 -tao_brgn_regularization_type l2prox -tao_brgn_regularizer_weight 1e-8 -tao_gatol 1.e-6
297c4762a1bSJed Brown 
298c4762a1bSJed Brown    test:
299c4762a1bSJed Brown       suffix: 3
300c4762a1bSJed Brown       localrunfiles: tomographyData_A_b_xGT
301c4762a1bSJed Brown       args: -tao_monitor -tao_max_it 1000 -tao_brgn_regularization_type user -tao_brgn_regularizer_weight 1e-8 -tao_gatol 1.e-6
302c4762a1bSJed Brown 
303c4762a1bSJed Brown TEST*/
304