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