xref: /petsc/src/tao/leastsquares/tutorials/chwirut1.c (revision cd1c4666a9c7786acffe2dfe971f76ff95fb06a8)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown    Include "petsctao.h" so that we can use TAO solvers.  Note that this
3c4762a1bSJed Brown    file automatically includes libraries such as:
4c4762a1bSJed Brown      petsc.h       - base PETSc routines   petscvec.h - vectors
5c4762a1bSJed Brown      petscsys.h    - sysem routines        petscmat.h - matrices
6c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
7c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
8c4762a1bSJed Brown 
9c4762a1bSJed Brown */
10c4762a1bSJed Brown 
11c4762a1bSJed Brown #include <petsctao.h>
12c4762a1bSJed Brown 
13c4762a1bSJed Brown /*
14c4762a1bSJed Brown Description:   These data are the result of a NIST study involving
15c4762a1bSJed Brown                ultrasonic calibration.  The response variable is
16c4762a1bSJed Brown                ultrasonic response, and the predictor variable is
17c4762a1bSJed Brown                metal distance.
18c4762a1bSJed Brown 
19c4762a1bSJed Brown Reference:     Chwirut, D., NIST (197?).
20c4762a1bSJed Brown                Ultrasonic Reference Block Study.
21c4762a1bSJed Brown */
22c4762a1bSJed Brown 
23c4762a1bSJed Brown static char help[]="Finds the nonlinear least-squares solution to the model \n\
24c4762a1bSJed Brown             y = exp[-b1*x]/(b2+b3*x)  +  e \n";
25c4762a1bSJed Brown 
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();
32c4762a1bSJed Brown    Routines: TaoSetInitialVector();
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 NOBSERVATIONS 214
41c4762a1bSJed Brown #define NPARAMETERS 3
42c4762a1bSJed Brown 
43c4762a1bSJed Brown /* User-defined application context */
44c4762a1bSJed Brown typedef struct {
45c4762a1bSJed Brown   /* Working space */
46c4762a1bSJed Brown   PetscReal t[NOBSERVATIONS];   /* array of independent variables of observation */
47c4762a1bSJed Brown   PetscReal y[NOBSERVATIONS];   /* array of dependent variables */
48c4762a1bSJed Brown   PetscReal j[NOBSERVATIONS][NPARAMETERS]; /* dense jacobian matrix array*/
49c4762a1bSJed Brown   PetscInt idm[NOBSERVATIONS];  /* Matrix indices for jacobian */
50c4762a1bSJed Brown   PetscInt idn[NPARAMETERS];
51c4762a1bSJed Brown } AppCtx;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown /* User provided Routines */
54c4762a1bSJed Brown PetscErrorCode InitializeData(AppCtx *user);
55c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec);
56c4762a1bSJed Brown PetscErrorCode EvaluateFunction(Tao, Vec, Vec, void *);
57c4762a1bSJed Brown PetscErrorCode EvaluateJacobian(Tao, Vec, Mat, Mat, void *);
58c4762a1bSJed Brown 
59c4762a1bSJed Brown /*--------------------------------------------------------------------*/
60c4762a1bSJed Brown int main(int argc,char **argv)
61c4762a1bSJed Brown {
62c4762a1bSJed Brown   PetscErrorCode ierr;           /* used to check for functions returning nonzeros */
63c4762a1bSJed Brown   Vec            x, f;               /* solution, function */
64c4762a1bSJed Brown   Mat            J;                  /* Jacobian matrix */
65c4762a1bSJed Brown   Tao            tao;                /* Tao solver context */
66c4762a1bSJed Brown   PetscInt       i;               /* iteration information */
67c4762a1bSJed Brown   PetscReal      hist[100],resid[100];
68c4762a1bSJed Brown   PetscInt       lits[100];
69c4762a1bSJed Brown   AppCtx         user;               /* user-defined work context */
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char *)0,help);if (ierr) return ierr;
72c4762a1bSJed Brown   /* Allocate vectors */
73c4762a1bSJed Brown   ierr = VecCreateSeq(MPI_COMM_SELF,NPARAMETERS,&x);CHKERRQ(ierr);
74c4762a1bSJed Brown   ierr = VecCreateSeq(MPI_COMM_SELF,NOBSERVATIONS,&f);CHKERRQ(ierr);
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* Create the Jacobian matrix. */
77c4762a1bSJed Brown   ierr = MatCreateSeqDense(MPI_COMM_SELF,NOBSERVATIONS,NPARAMETERS,NULL,&J);CHKERRQ(ierr);
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   for (i=0;i<NOBSERVATIONS;i++) user.idm[i] = i;
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   for (i=0;i<NPARAMETERS;i++) user.idn[i] = i;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
84c4762a1bSJed Brown   ierr = TaoCreate(PETSC_COMM_SELF,&tao);CHKERRQ(ierr);
85c4762a1bSJed Brown   ierr = TaoSetType(tao,TAOPOUNDERS);CHKERRQ(ierr);
86c4762a1bSJed Brown 
87c4762a1bSJed Brown  /* Set the function and Jacobian routines. */
88c4762a1bSJed Brown   ierr = InitializeData(&user);CHKERRQ(ierr);
89c4762a1bSJed Brown   ierr = FormStartingPoint(x);CHKERRQ(ierr);
90c4762a1bSJed Brown   ierr = TaoSetInitialVector(tao,x);CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = TaoSetResidualRoutine(tao,f,EvaluateFunction,(void*)&user);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = TaoSetJacobianResidualRoutine(tao, J, J, EvaluateJacobian, (void*)&user);CHKERRQ(ierr);
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /* Check for any TAO command line arguments */
95c4762a1bSJed Brown   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   ierr = TaoSetConvergenceHistory(tao,hist,resid,0,lits,100,PETSC_TRUE);CHKERRQ(ierr);
98c4762a1bSJed Brown   /* Perform the Solve */
99c4762a1bSJed Brown   ierr = TaoSolve(tao);CHKERRQ(ierr);
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /* View the vector; then destroy it.  */
102c4762a1bSJed Brown   ierr = VecView(x,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   /* Free TAO data structures */
105c4762a1bSJed Brown   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
106c4762a1bSJed Brown 
107c4762a1bSJed Brown    /* Free PETSc data structures */
108c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
109c4762a1bSJed Brown   ierr = VecDestroy(&f);CHKERRQ(ierr);
110c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   ierr = PetscFinalize();
113c4762a1bSJed Brown   return ierr;
114c4762a1bSJed Brown }
115c4762a1bSJed Brown 
116c4762a1bSJed Brown /*--------------------------------------------------------------------*/
117c4762a1bSJed Brown PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr)
118c4762a1bSJed Brown {
119c4762a1bSJed Brown   AppCtx         *user = (AppCtx *)ptr;
120c4762a1bSJed Brown   PetscInt       i;
121c4762a1bSJed Brown   const PetscReal *x;
122c4762a1bSJed Brown   PetscReal      *y=user->y,*f,*t=user->t;
123c4762a1bSJed Brown   PetscErrorCode ierr;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   PetscFunctionBegin;
126c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
127c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   for (i=0;i<NOBSERVATIONS;i++) {
130c4762a1bSJed Brown     f[i] = y[i] - PetscExpScalar(-x[0]*t[i])/(x[1] + x[2]*t[i]);
131c4762a1bSJed Brown   }
132c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
133c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
134c4762a1bSJed Brown   PetscLogFlops(6*NOBSERVATIONS);
135c4762a1bSJed Brown   PetscFunctionReturn(0);
136c4762a1bSJed Brown }
137c4762a1bSJed Brown 
138c4762a1bSJed Brown /*------------------------------------------------------------*/
139c4762a1bSJed Brown /* J[i][j] = df[i]/dt[j] */
140c4762a1bSJed Brown PetscErrorCode EvaluateJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr)
141c4762a1bSJed Brown {
142c4762a1bSJed Brown   AppCtx         *user = (AppCtx *)ptr;
143c4762a1bSJed Brown   PetscInt       i;
144c4762a1bSJed Brown   const PetscReal *x;
145c4762a1bSJed Brown   PetscReal      *t=user->t;
146c4762a1bSJed Brown   PetscReal      base;
147c4762a1bSJed Brown   PetscErrorCode ierr;
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   PetscFunctionBegin;
150c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
151c4762a1bSJed Brown   for (i=0;i<NOBSERVATIONS;i++) {
152c4762a1bSJed Brown     base = PetscExpScalar(-x[0]*t[i])/(x[1] + x[2]*t[i]);
153c4762a1bSJed Brown 
154c4762a1bSJed Brown     user->j[i][0] = t[i]*base;
155c4762a1bSJed Brown     user->j[i][1] = base/(x[1] + x[2]*t[i]);
156c4762a1bSJed Brown     user->j[i][2] = base*t[i]/(x[1] + x[2]*t[i]);
157c4762a1bSJed Brown   }
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /* Assemble the matrix */
160c4762a1bSJed Brown   ierr = MatSetValues(J,NOBSERVATIONS,user->idm, NPARAMETERS, user->idn,(PetscReal *)user->j,INSERT_VALUES);CHKERRQ(ierr);
161c4762a1bSJed Brown   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
162c4762a1bSJed Brown   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
165c4762a1bSJed Brown   PetscLogFlops(NOBSERVATIONS * 13);
166c4762a1bSJed Brown   PetscFunctionReturn(0);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown /* ------------------------------------------------------------ */
170c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec X)
171c4762a1bSJed Brown {
172c4762a1bSJed Brown   PetscReal      *x;
173c4762a1bSJed Brown   PetscErrorCode ierr;
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   PetscFunctionBegin;
176c4762a1bSJed Brown   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
177c4762a1bSJed Brown   x[0] = 0.15;
178c4762a1bSJed Brown   x[1] = 0.008;
179c4762a1bSJed Brown   x[2] = 0.010;
180c4762a1bSJed Brown   VecRestoreArray(X,&x);CHKERRQ(ierr);
181c4762a1bSJed Brown   PetscFunctionReturn(0);
182c4762a1bSJed Brown }
183c4762a1bSJed Brown 
184c4762a1bSJed Brown /* ---------------------------------------------------------------------- */
185c4762a1bSJed Brown PetscErrorCode InitializeData(AppCtx *user)
186c4762a1bSJed Brown {
187c4762a1bSJed Brown   PetscReal *t=user->t,*y=user->y;
188c4762a1bSJed Brown   PetscInt  i=0;
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   PetscFunctionBegin;
191c4762a1bSJed Brown   y[i] =   92.9000;   t[i++] =  0.5000;
192c4762a1bSJed Brown   y[i] =    78.7000;  t[i++] =   0.6250;
193c4762a1bSJed Brown   y[i] =    64.2000;  t[i++] =   0.7500;
194c4762a1bSJed Brown   y[i] =    64.9000;  t[i++] =   0.8750;
195c4762a1bSJed Brown   y[i] =    57.1000;  t[i++] =   1.0000;
196c4762a1bSJed Brown   y[i] =    43.3000;  t[i++] =   1.2500;
197c4762a1bSJed Brown   y[i] =    31.1000;   t[i++] =  1.7500;
198c4762a1bSJed Brown   y[i] =    23.6000;   t[i++] =  2.2500;
199c4762a1bSJed Brown   y[i] =    31.0500;   t[i++] =  1.7500;
200c4762a1bSJed Brown   y[i] =    23.7750;   t[i++] =  2.2500;
201c4762a1bSJed Brown   y[i] =    17.7375;   t[i++] =  2.7500;
202c4762a1bSJed Brown   y[i] =    13.8000;   t[i++] =  3.2500;
203c4762a1bSJed Brown   y[i] =    11.5875;   t[i++] =  3.7500;
204c4762a1bSJed Brown   y[i] =     9.4125;   t[i++] =  4.2500;
205c4762a1bSJed Brown   y[i] =     7.7250;   t[i++] =  4.7500;
206c4762a1bSJed Brown   y[i] =     7.3500;   t[i++] =  5.2500;
207c4762a1bSJed Brown   y[i] =     8.0250;   t[i++] =  5.7500;
208c4762a1bSJed Brown   y[i] =    90.6000;   t[i++] =  0.5000;
209c4762a1bSJed Brown   y[i] =    76.9000;   t[i++] =  0.6250;
210c4762a1bSJed Brown   y[i] =    71.6000;   t[i++] = 0.7500;
211c4762a1bSJed Brown   y[i] =    63.6000;   t[i++] =  0.8750;
212c4762a1bSJed Brown   y[i] =    54.0000;   t[i++] =  1.0000;
213c4762a1bSJed Brown   y[i] =    39.2000;   t[i++] =  1.2500;
214c4762a1bSJed Brown   y[i] =    29.3000;   t[i++] = 1.7500;
215c4762a1bSJed Brown   y[i] =    21.4000;   t[i++] =  2.2500;
216c4762a1bSJed Brown   y[i] =    29.1750;   t[i++] =  1.7500;
217c4762a1bSJed Brown   y[i] =    22.1250;   t[i++] =  2.2500;
218c4762a1bSJed Brown   y[i] =    17.5125;   t[i++] =  2.7500;
219c4762a1bSJed Brown   y[i] =    14.2500;   t[i++] =  3.2500;
220c4762a1bSJed Brown   y[i] =     9.4500;   t[i++] =  3.7500;
221c4762a1bSJed Brown   y[i] =     9.1500;   t[i++] =  4.2500;
222c4762a1bSJed Brown   y[i] =     7.9125;   t[i++] =  4.7500;
223c4762a1bSJed Brown   y[i] =     8.4750;   t[i++] =  5.2500;
224c4762a1bSJed Brown   y[i] =     6.1125;   t[i++] =  5.7500;
225c4762a1bSJed Brown   y[i] =    80.0000;   t[i++] =  0.5000;
226c4762a1bSJed Brown   y[i] =    79.0000;   t[i++] =  0.6250;
227c4762a1bSJed Brown   y[i] =    63.8000;   t[i++] =  0.7500;
228c4762a1bSJed Brown   y[i] =    57.2000;   t[i++] =  0.8750;
229c4762a1bSJed Brown   y[i] =    53.2000;   t[i++] =  1.0000;
230c4762a1bSJed Brown   y[i] =   42.5000;   t[i++] =  1.2500;
231c4762a1bSJed Brown   y[i] =   26.8000;   t[i++] =  1.7500;
232c4762a1bSJed Brown   y[i] =    20.4000;   t[i++] =  2.2500;
233c4762a1bSJed Brown   y[i] =    26.8500;  t[i++] =   1.7500;
234c4762a1bSJed Brown   y[i] =    21.0000;  t[i++] =   2.2500;
235c4762a1bSJed Brown   y[i] =    16.4625;  t[i++] =   2.7500;
236c4762a1bSJed Brown   y[i] =    12.5250;  t[i++] =   3.2500;
237c4762a1bSJed Brown   y[i] =    10.5375;  t[i++] =   3.7500;
238c4762a1bSJed Brown   y[i] =     8.5875;  t[i++] =   4.2500;
239c4762a1bSJed Brown   y[i] =     7.1250;  t[i++] =   4.7500;
240c4762a1bSJed Brown   y[i] =     6.1125;  t[i++] =   5.2500;
241c4762a1bSJed Brown   y[i] =     5.9625;  t[i++] =   5.7500;
242c4762a1bSJed Brown   y[i] =    74.1000;  t[i++] =   0.5000;
243c4762a1bSJed Brown   y[i] =    67.3000;  t[i++] =   0.6250;
244c4762a1bSJed Brown   y[i] =    60.8000;  t[i++] =   0.7500;
245c4762a1bSJed Brown   y[i] =    55.5000;  t[i++] =   0.8750;
246c4762a1bSJed Brown   y[i] =    50.3000;  t[i++] =   1.0000;
247c4762a1bSJed Brown   y[i] =    41.0000;  t[i++] =   1.2500;
248c4762a1bSJed Brown   y[i] =    29.4000;  t[i++] =   1.7500;
249c4762a1bSJed Brown   y[i] =    20.4000;  t[i++] =   2.2500;
250c4762a1bSJed Brown   y[i] =    29.3625;  t[i++] =   1.7500;
251c4762a1bSJed Brown   y[i] =    21.1500;  t[i++] =   2.2500;
252c4762a1bSJed Brown   y[i] =    16.7625;  t[i++] =   2.7500;
253c4762a1bSJed Brown   y[i] =    13.2000;  t[i++] =   3.2500;
254c4762a1bSJed Brown   y[i] =    10.8750;  t[i++] =   3.7500;
255c4762a1bSJed Brown   y[i] =     8.1750;  t[i++] =   4.2500;
256c4762a1bSJed Brown   y[i] =     7.3500;  t[i++] =   4.7500;
257c4762a1bSJed Brown   y[i] =     5.9625;  t[i++] =  5.2500;
258c4762a1bSJed Brown   y[i] =     5.6250;  t[i++] =   5.7500;
259c4762a1bSJed Brown   y[i] =    81.5000;  t[i++] =    .5000;
260c4762a1bSJed Brown   y[i] =    62.4000;  t[i++] =    .7500;
261c4762a1bSJed Brown   y[i] =    32.5000;  t[i++] =   1.5000;
262c4762a1bSJed Brown   y[i] =    12.4100;  t[i++] =   3.0000;
263c4762a1bSJed Brown   y[i] =    13.1200;  t[i++] =   3.0000;
264c4762a1bSJed Brown   y[i] =    15.5600;  t[i++] =   3.0000;
265c4762a1bSJed Brown   y[i] =     5.6300;  t[i++] =   6.0000;
266c4762a1bSJed Brown   y[i] =    78.0000;   t[i++] =   .5000;
267c4762a1bSJed Brown   y[i] =    59.9000;  t[i++] =    .7500;
268c4762a1bSJed Brown   y[i] =    33.2000;  t[i++] =   1.5000;
269c4762a1bSJed Brown   y[i] =    13.8400;  t[i++] =   3.0000;
270c4762a1bSJed Brown   y[i] =    12.7500;  t[i++] =   3.0000;
271c4762a1bSJed Brown   y[i] =    14.6200;  t[i++] =   3.0000;
272c4762a1bSJed Brown   y[i] =     3.9400;  t[i++] =   6.0000;
273c4762a1bSJed Brown   y[i] =    76.8000;  t[i++] =    .5000;
274c4762a1bSJed Brown   y[i] =    61.0000;  t[i++] =    .7500;
275c4762a1bSJed Brown   y[i] =    32.9000;  t[i++] =   1.5000;
276c4762a1bSJed Brown   y[i] =   13.8700;   t[i++] = 3.0000;
277c4762a1bSJed Brown   y[i] =    11.8100;  t[i++] =   3.0000;
278c4762a1bSJed Brown   y[i] =    13.3100;  t[i++] =   3.0000;
279c4762a1bSJed Brown   y[i] =     5.4400;  t[i++] =   6.0000;
280c4762a1bSJed Brown   y[i] =    78.0000;  t[i++] =    .5000;
281c4762a1bSJed Brown   y[i] =    63.5000;  t[i++] =    .7500;
282c4762a1bSJed Brown   y[i] =    33.8000;  t[i++] =   1.5000;
283c4762a1bSJed Brown   y[i] =    12.5600;  t[i++] =   3.0000;
284c4762a1bSJed Brown   y[i] =     5.6300;  t[i++] =   6.0000;
285c4762a1bSJed Brown   y[i] =    12.7500;  t[i++] =   3.0000;
286c4762a1bSJed Brown   y[i] =    13.1200;  t[i++] =   3.0000;
287c4762a1bSJed Brown   y[i] =     5.4400;  t[i++] =   6.0000;
288c4762a1bSJed Brown   y[i] =    76.8000;  t[i++] =    .5000;
289c4762a1bSJed Brown   y[i] =    60.0000;  t[i++] =    .7500;
290c4762a1bSJed Brown   y[i] =    47.8000;  t[i++] =   1.0000;
291c4762a1bSJed Brown   y[i] =    32.0000;  t[i++] =   1.5000;
292c4762a1bSJed Brown   y[i] =    22.2000;  t[i++] =   2.0000;
293c4762a1bSJed Brown   y[i] =    22.5700;  t[i++] =   2.0000;
294c4762a1bSJed Brown   y[i] =    18.8200;  t[i++] =   2.5000;
295c4762a1bSJed Brown   y[i] =    13.9500;  t[i++] =   3.0000;
296c4762a1bSJed Brown   y[i] =    11.2500;  t[i++] =   4.0000;
297c4762a1bSJed Brown   y[i] =     9.0000;  t[i++] =   5.0000;
298c4762a1bSJed Brown   y[i] =     6.6700;  t[i++] =   6.0000;
299c4762a1bSJed Brown   y[i] =    75.8000;  t[i++] =    .5000;
300c4762a1bSJed Brown   y[i] =    62.0000;  t[i++] =    .7500;
301c4762a1bSJed Brown   y[i] =    48.8000;  t[i++] =   1.0000;
302c4762a1bSJed Brown   y[i] =    35.2000;  t[i++] =   1.5000;
303c4762a1bSJed Brown   y[i] =    20.0000;  t[i++] =   2.0000;
304c4762a1bSJed Brown   y[i] =    20.3200;  t[i++] =   2.0000;
305c4762a1bSJed Brown   y[i] =    19.3100;  t[i++] =   2.5000;
306c4762a1bSJed Brown   y[i] =    12.7500;  t[i++] =   3.0000;
307c4762a1bSJed Brown   y[i] =    10.4200;  t[i++] =   4.0000;
308c4762a1bSJed Brown   y[i] =     7.3100;  t[i++] =   5.0000;
309c4762a1bSJed Brown   y[i] =     7.4200;  t[i++] =   6.0000;
310c4762a1bSJed Brown   y[i] =    70.5000;  t[i++] =    .5000;
311c4762a1bSJed Brown   y[i] =    59.5000;  t[i++] =    .7500;
312c4762a1bSJed Brown   y[i] =    48.5000;  t[i++] =   1.0000;
313c4762a1bSJed Brown   y[i] =    35.8000;  t[i++] =   1.5000;
314c4762a1bSJed Brown   y[i] =    21.0000;  t[i++] =   2.0000;
315c4762a1bSJed Brown   y[i] =    21.6700;  t[i++] =   2.0000;
316c4762a1bSJed Brown   y[i] =    21.0000;  t[i++] =   2.5000;
317c4762a1bSJed Brown   y[i] =    15.6400;  t[i++] =   3.0000;
318c4762a1bSJed Brown   y[i] =     8.1700;  t[i++] =   4.0000;
319c4762a1bSJed Brown   y[i] =     8.5500;  t[i++] =   5.0000;
320c4762a1bSJed Brown   y[i] =    10.1200;  t[i++] =   6.0000;
321c4762a1bSJed Brown   y[i] =    78.0000;  t[i++] =    .5000;
322c4762a1bSJed Brown   y[i] =    66.0000;  t[i++] =    .6250;
323c4762a1bSJed Brown   y[i] =    62.0000;  t[i++] =    .7500;
324c4762a1bSJed Brown   y[i] =    58.0000;  t[i++] =    .8750;
325c4762a1bSJed Brown   y[i] =    47.7000;  t[i++] =   1.0000;
326c4762a1bSJed Brown   y[i] =    37.8000;  t[i++] =   1.2500;
327c4762a1bSJed Brown   y[i] =    20.2000;  t[i++] =   2.2500;
328c4762a1bSJed Brown   y[i] =    21.0700;  t[i++] =   2.2500;
329c4762a1bSJed Brown   y[i] =    13.8700;  t[i++] =   2.7500;
330c4762a1bSJed Brown   y[i] =     9.6700;  t[i++] =   3.2500;
331c4762a1bSJed Brown   y[i] =     7.7600;  t[i++] =   3.7500;
332c4762a1bSJed Brown   y[i] =    5.4400;   t[i++] =  4.2500;
333c4762a1bSJed Brown   y[i] =    4.8700;   t[i++] =  4.7500;
334c4762a1bSJed Brown   y[i] =     4.0100;  t[i++] =   5.2500;
335c4762a1bSJed Brown   y[i] =     3.7500;  t[i++] =   5.7500;
336c4762a1bSJed Brown   y[i] =    24.1900;  t[i++] =   3.0000;
337c4762a1bSJed Brown   y[i] =    25.7600;  t[i++] =   3.0000;
338c4762a1bSJed Brown   y[i] =    18.0700;  t[i++] =   3.0000;
339c4762a1bSJed Brown   y[i] =    11.8100;  t[i++] =   3.0000;
340c4762a1bSJed Brown   y[i] =    12.0700;  t[i++] =   3.0000;
341c4762a1bSJed Brown   y[i] =    16.1200;  t[i++] =   3.0000;
342c4762a1bSJed Brown   y[i] =    70.8000;  t[i++] =    .5000;
343c4762a1bSJed Brown   y[i] =    54.7000;  t[i++] =    .7500;
344c4762a1bSJed Brown   y[i] =    48.0000;  t[i++] =   1.0000;
345c4762a1bSJed Brown   y[i] =    39.8000;  t[i++] =   1.5000;
346c4762a1bSJed Brown   y[i] =    29.8000;  t[i++] =   2.0000;
347c4762a1bSJed Brown   y[i] =    23.7000;  t[i++] =   2.5000;
348c4762a1bSJed Brown   y[i] =    29.6200;  t[i++] =   2.0000;
349c4762a1bSJed Brown   y[i] =    23.8100;  t[i++] =   2.5000;
350c4762a1bSJed Brown   y[i] =    17.7000;  t[i++] =   3.0000;
351c4762a1bSJed Brown   y[i] =    11.5500;  t[i++] =   4.0000;
352c4762a1bSJed Brown   y[i] =    12.0700;  t[i++] =   5.0000;
353c4762a1bSJed Brown   y[i] =     8.7400;  t[i++] =   6.0000;
354c4762a1bSJed Brown   y[i] =    80.7000;  t[i++] =    .5000;
355c4762a1bSJed Brown   y[i] =    61.3000;  t[i++] =    .7500;
356c4762a1bSJed Brown   y[i] =    47.5000;  t[i++] =   1.0000;
357c4762a1bSJed Brown    y[i] =   29.0000;  t[i++] =   1.5000;
358c4762a1bSJed Brown    y[i] =   24.0000;  t[i++] =   2.0000;
359c4762a1bSJed Brown   y[i] =    17.7000;  t[i++] =   2.5000;
360c4762a1bSJed Brown   y[i] =    24.5600;  t[i++] =   2.0000;
361c4762a1bSJed Brown   y[i] =    18.6700;  t[i++] =   2.5000;
362c4762a1bSJed Brown    y[i] =   16.2400;  t[i++] =   3.0000;
363c4762a1bSJed Brown   y[i] =     8.7400;  t[i++] =   4.0000;
364c4762a1bSJed Brown   y[i] =     7.8700;  t[i++] =   5.0000;
365c4762a1bSJed Brown   y[i] =     8.5100;  t[i++] =   6.0000;
366c4762a1bSJed Brown   y[i] =    66.7000;  t[i++] =    .5000;
367c4762a1bSJed Brown   y[i] =    59.2000;  t[i++] =    .7500;
368c4762a1bSJed Brown   y[i] =    40.8000;  t[i++] =   1.0000;
369c4762a1bSJed Brown   y[i] =    30.7000;  t[i++] =   1.5000;
370c4762a1bSJed Brown   y[i] =    25.7000;  t[i++] =   2.0000;
371c4762a1bSJed Brown   y[i] =    16.3000;  t[i++] =   2.5000;
372c4762a1bSJed Brown   y[i] =    25.9900;  t[i++] =   2.0000;
373c4762a1bSJed Brown   y[i] =    16.9500;  t[i++] =   2.5000;
374c4762a1bSJed Brown   y[i] =    13.3500;  t[i++] =   3.0000;
375c4762a1bSJed Brown   y[i] =     8.6200;  t[i++] =   4.0000;
376c4762a1bSJed Brown   y[i] =     7.2000;  t[i++] =   5.0000;
377c4762a1bSJed Brown   y[i] =     6.6400;  t[i++] =   6.0000;
378c4762a1bSJed Brown   y[i] =    13.6900;  t[i++] =   3.0000;
379c4762a1bSJed Brown   y[i] =    81.0000;  t[i++] =    .5000;
380c4762a1bSJed Brown   y[i] =    64.5000;  t[i++] =    .7500;
381c4762a1bSJed Brown   y[i] =    35.5000;  t[i++] =   1.5000;
382c4762a1bSJed Brown    y[i] =   13.3100;  t[i++] =   3.0000;
383c4762a1bSJed Brown   y[i] =     4.8700;  t[i++] =   6.0000;
384c4762a1bSJed Brown   y[i] =    12.9400;  t[i++] =   3.0000;
385c4762a1bSJed Brown   y[i] =     5.0600;  t[i++] =   6.0000;
386c4762a1bSJed Brown   y[i] =    15.1900;  t[i++] =   3.0000;
387c4762a1bSJed Brown   y[i] =    14.6200;  t[i++] =   3.0000;
388c4762a1bSJed Brown   y[i] =    15.6400;  t[i++] =   3.0000;
389c4762a1bSJed Brown   y[i] =    25.5000;  t[i++] =   1.7500;
390c4762a1bSJed Brown   y[i] =    25.9500;  t[i++] =   1.7500;
391c4762a1bSJed Brown   y[i] =    81.7000;  t[i++] =    .5000;
392c4762a1bSJed Brown   y[i] =    61.6000;  t[i++] =    .7500;
393c4762a1bSJed Brown   y[i] =    29.8000;  t[i++] =   1.7500;
394c4762a1bSJed Brown   y[i] =    29.8100;  t[i++] =   1.7500;
395c4762a1bSJed Brown   y[i] =    17.1700;  t[i++] =   2.7500;
396c4762a1bSJed Brown   y[i] =    10.3900;  t[i++] =   3.7500;
397c4762a1bSJed Brown   y[i] =    28.4000;  t[i++] =   1.7500;
398c4762a1bSJed Brown   y[i] =    28.6900;  t[i++] =   1.7500;
399c4762a1bSJed Brown   y[i] =    81.3000;  t[i++] =    .5000;
400c4762a1bSJed Brown   y[i] =    60.9000;  t[i++] =    .7500;
401c4762a1bSJed Brown   y[i] =    16.6500;  t[i++] =   2.7500;
402c4762a1bSJed Brown   y[i] =    10.0500;  t[i++] =   3.7500;
403c4762a1bSJed Brown   y[i] =    28.9000;  t[i++] =   1.7500;
404c4762a1bSJed Brown   y[i] =    28.9500;  t[i++] =   1.7500;
405c4762a1bSJed Brown   PetscFunctionReturn(0);
406c4762a1bSJed Brown }
407c4762a1bSJed Brown 
408c4762a1bSJed Brown /*TEST
409c4762a1bSJed Brown 
410c4762a1bSJed Brown    build:
411c4762a1bSJed Brown       requires: !complex !single
412c4762a1bSJed Brown 
413c4762a1bSJed Brown    test:
414c4762a1bSJed Brown       args: -tao_smonitor -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-5
415c4762a1bSJed Brown 
416c4762a1bSJed Brown    test:
417c4762a1bSJed Brown       suffix: 2
418c4762a1bSJed Brown       args: -tao_smonitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type l2prox -tao_brgn_regularizer_weight 1e-4 -tao_gatol 1.e-5
419c4762a1bSJed Brown 
420c4762a1bSJed Brown    test:
421c4762a1bSJed Brown       suffix: 3
422c4762a1bSJed Brown       args: -tao_smonitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type l1dict -tao_brgn_regularizer_weight 1e-4 -tao_brgn_l1_smooth_epsilon 1e-6 -tao_gatol 1.e-5
423c4762a1bSJed Brown 
424*cd1c4666STristan Konolige    test:
425*cd1c4666STristan Konolige       suffix: 4
426*cd1c4666STristan Konolige       args: -tao_smonitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type lm -tao_gatol 1.e-5 -tao_brgn_subsolver_tao_type bnls
427*cd1c4666STristan Konolige 
428c4762a1bSJed Brown TEST*/
429