xref: /petsc/src/tao/leastsquares/tutorials/chwirut1.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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
5a5b23f4aSJose E. Roman      petscsys.h    - system 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 #define NOBSERVATIONS 214
27c4762a1bSJed Brown #define NPARAMETERS   3
28c4762a1bSJed Brown 
29c4762a1bSJed Brown /* User-defined application context */
30c4762a1bSJed Brown typedef struct {
31c4762a1bSJed Brown   /* Working space */
32c4762a1bSJed Brown   PetscReal t[NOBSERVATIONS];              /* array of independent variables of observation */
33c4762a1bSJed Brown   PetscReal y[NOBSERVATIONS];              /* array of dependent variables */
34c4762a1bSJed Brown   PetscReal j[NOBSERVATIONS][NPARAMETERS]; /* dense jacobian matrix array*/
35c4762a1bSJed Brown   PetscInt  idm[NOBSERVATIONS];            /* Matrix indices for jacobian */
36c4762a1bSJed Brown   PetscInt  idn[NPARAMETERS];
37c4762a1bSJed Brown } AppCtx;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown /* User provided Routines */
40c4762a1bSJed Brown PetscErrorCode InitializeData(AppCtx *user);
41c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec);
42c4762a1bSJed Brown PetscErrorCode EvaluateFunction(Tao, Vec, Vec, void *);
43c4762a1bSJed Brown PetscErrorCode EvaluateJacobian(Tao, Vec, Mat, Mat, void *);
44c4762a1bSJed Brown 
45c4762a1bSJed Brown /*--------------------------------------------------------------------*/
46*9371c9d4SSatish Balay int main(int argc, char **argv) {
47c4762a1bSJed Brown   Vec       x, f; /* solution, function */
48c4762a1bSJed Brown   Mat       J;    /* Jacobian matrix */
49c4762a1bSJed Brown   Tao       tao;  /* Tao solver context */
50c4762a1bSJed Brown   PetscInt  i;    /* iteration information */
51c4762a1bSJed Brown   PetscReal hist[100], resid[100];
52c4762a1bSJed Brown   PetscInt  lits[100];
53c4762a1bSJed Brown   AppCtx    user; /* user-defined work context */
54c4762a1bSJed Brown 
55327415f7SBarry Smith   PetscFunctionBeginUser;
569566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
57c4762a1bSJed Brown   /* Allocate vectors */
589566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(MPI_COMM_SELF, NPARAMETERS, &x));
599566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(MPI_COMM_SELF, NOBSERVATIONS, &f));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Create the Jacobian matrix. */
629566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(MPI_COMM_SELF, NOBSERVATIONS, NPARAMETERS, NULL, &J));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   for (i = 0; i < NOBSERVATIONS; i++) user.idm[i] = i;
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   for (i = 0; i < NPARAMETERS; i++) user.idn[i] = i;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
699566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_SELF, &tao));
709566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOPOUNDERS));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* Set the function and Jacobian routines. */
739566063dSJacob Faibussowitsch   PetscCall(InitializeData(&user));
749566063dSJacob Faibussowitsch   PetscCall(FormStartingPoint(x));
759566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, x));
769566063dSJacob Faibussowitsch   PetscCall(TaoSetResidualRoutine(tao, f, EvaluateFunction, (void *)&user));
779566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianResidualRoutine(tao, J, J, EvaluateJacobian, (void *)&user));
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* Check for any TAO command line arguments */
809566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
81c4762a1bSJed Brown 
829566063dSJacob Faibussowitsch   PetscCall(TaoSetConvergenceHistory(tao, hist, resid, 0, lits, 100, PETSC_TRUE));
83c4762a1bSJed Brown   /* Perform the Solve */
849566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* View the vector; then destroy it.  */
879566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_SELF));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* Free TAO data structures */
909566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   /* Free PETSc data structures */
939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&f));
959566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
96c4762a1bSJed Brown 
979566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
98b122ec5aSJacob Faibussowitsch   return 0;
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown /*--------------------------------------------------------------------*/
102*9371c9d4SSatish Balay PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr) {
103c4762a1bSJed Brown   AppCtx          *user = (AppCtx *)ptr;
104c4762a1bSJed Brown   PetscInt         i;
105c4762a1bSJed Brown   const PetscReal *x;
106c4762a1bSJed Brown   PetscReal       *y = user->y, *f, *t = user->t;
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   PetscFunctionBegin;
1099566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1109566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
111c4762a1bSJed Brown 
112*9371c9d4SSatish Balay   for (i = 0; i < NOBSERVATIONS; i++) { f[i] = y[i] - PetscExpScalar(-x[0] * t[i]) / (x[1] + x[2] * t[i]); }
1139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1149566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
115c4762a1bSJed Brown   PetscLogFlops(6 * NOBSERVATIONS);
116c4762a1bSJed Brown   PetscFunctionReturn(0);
117c4762a1bSJed Brown }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown /*------------------------------------------------------------*/
120c4762a1bSJed Brown /* J[i][j] = df[i]/dt[j] */
121*9371c9d4SSatish Balay PetscErrorCode EvaluateJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr) {
122c4762a1bSJed Brown   AppCtx          *user = (AppCtx *)ptr;
123c4762a1bSJed Brown   PetscInt         i;
124c4762a1bSJed Brown   const PetscReal *x;
125c4762a1bSJed Brown   PetscReal       *t = user->t;
126c4762a1bSJed Brown   PetscReal        base;
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   PetscFunctionBegin;
1299566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
130c4762a1bSJed Brown   for (i = 0; i < NOBSERVATIONS; i++) {
131c4762a1bSJed Brown     base = PetscExpScalar(-x[0] * t[i]) / (x[1] + x[2] * t[i]);
132c4762a1bSJed Brown 
133c4762a1bSJed Brown     user->j[i][0] = t[i] * base;
134c4762a1bSJed Brown     user->j[i][1] = base / (x[1] + x[2] * t[i]);
135c4762a1bSJed Brown     user->j[i][2] = base * t[i] / (x[1] + x[2] * t[i]);
136c4762a1bSJed Brown   }
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   /* Assemble the matrix */
1399566063dSJacob Faibussowitsch   PetscCall(MatSetValues(J, NOBSERVATIONS, user->idm, NPARAMETERS, user->idn, (PetscReal *)user->j, INSERT_VALUES));
1409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
142c4762a1bSJed Brown 
1439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
144c4762a1bSJed Brown   PetscLogFlops(NOBSERVATIONS * 13);
145c4762a1bSJed Brown   PetscFunctionReturn(0);
146c4762a1bSJed Brown }
147c4762a1bSJed Brown 
148c4762a1bSJed Brown /* ------------------------------------------------------------ */
149*9371c9d4SSatish Balay PetscErrorCode FormStartingPoint(Vec X) {
150c4762a1bSJed Brown   PetscReal *x;
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   PetscFunctionBegin;
1539566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
154c4762a1bSJed Brown   x[0] = 0.15;
155c4762a1bSJed Brown   x[1] = 0.008;
156c4762a1bSJed Brown   x[2] = 0.010;
1579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
158c4762a1bSJed Brown   PetscFunctionReturn(0);
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown /* ---------------------------------------------------------------------- */
162*9371c9d4SSatish Balay PetscErrorCode InitializeData(AppCtx *user) {
163c4762a1bSJed Brown   PetscReal *t = user->t, *y = user->y;
164c4762a1bSJed Brown   PetscInt   i = 0;
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   PetscFunctionBegin;
167*9371c9d4SSatish Balay   y[i]   = 92.9000;
168*9371c9d4SSatish Balay   t[i++] = 0.5000;
169*9371c9d4SSatish Balay   y[i]   = 78.7000;
170*9371c9d4SSatish Balay   t[i++] = 0.6250;
171*9371c9d4SSatish Balay   y[i]   = 64.2000;
172*9371c9d4SSatish Balay   t[i++] = 0.7500;
173*9371c9d4SSatish Balay   y[i]   = 64.9000;
174*9371c9d4SSatish Balay   t[i++] = 0.8750;
175*9371c9d4SSatish Balay   y[i]   = 57.1000;
176*9371c9d4SSatish Balay   t[i++] = 1.0000;
177*9371c9d4SSatish Balay   y[i]   = 43.3000;
178*9371c9d4SSatish Balay   t[i++] = 1.2500;
179*9371c9d4SSatish Balay   y[i]   = 31.1000;
180*9371c9d4SSatish Balay   t[i++] = 1.7500;
181*9371c9d4SSatish Balay   y[i]   = 23.6000;
182*9371c9d4SSatish Balay   t[i++] = 2.2500;
183*9371c9d4SSatish Balay   y[i]   = 31.0500;
184*9371c9d4SSatish Balay   t[i++] = 1.7500;
185*9371c9d4SSatish Balay   y[i]   = 23.7750;
186*9371c9d4SSatish Balay   t[i++] = 2.2500;
187*9371c9d4SSatish Balay   y[i]   = 17.7375;
188*9371c9d4SSatish Balay   t[i++] = 2.7500;
189*9371c9d4SSatish Balay   y[i]   = 13.8000;
190*9371c9d4SSatish Balay   t[i++] = 3.2500;
191*9371c9d4SSatish Balay   y[i]   = 11.5875;
192*9371c9d4SSatish Balay   t[i++] = 3.7500;
193*9371c9d4SSatish Balay   y[i]   = 9.4125;
194*9371c9d4SSatish Balay   t[i++] = 4.2500;
195*9371c9d4SSatish Balay   y[i]   = 7.7250;
196*9371c9d4SSatish Balay   t[i++] = 4.7500;
197*9371c9d4SSatish Balay   y[i]   = 7.3500;
198*9371c9d4SSatish Balay   t[i++] = 5.2500;
199*9371c9d4SSatish Balay   y[i]   = 8.0250;
200*9371c9d4SSatish Balay   t[i++] = 5.7500;
201*9371c9d4SSatish Balay   y[i]   = 90.6000;
202*9371c9d4SSatish Balay   t[i++] = 0.5000;
203*9371c9d4SSatish Balay   y[i]   = 76.9000;
204*9371c9d4SSatish Balay   t[i++] = 0.6250;
205*9371c9d4SSatish Balay   y[i]   = 71.6000;
206*9371c9d4SSatish Balay   t[i++] = 0.7500;
207*9371c9d4SSatish Balay   y[i]   = 63.6000;
208*9371c9d4SSatish Balay   t[i++] = 0.8750;
209*9371c9d4SSatish Balay   y[i]   = 54.0000;
210*9371c9d4SSatish Balay   t[i++] = 1.0000;
211*9371c9d4SSatish Balay   y[i]   = 39.2000;
212*9371c9d4SSatish Balay   t[i++] = 1.2500;
213*9371c9d4SSatish Balay   y[i]   = 29.3000;
214*9371c9d4SSatish Balay   t[i++] = 1.7500;
215*9371c9d4SSatish Balay   y[i]   = 21.4000;
216*9371c9d4SSatish Balay   t[i++] = 2.2500;
217*9371c9d4SSatish Balay   y[i]   = 29.1750;
218*9371c9d4SSatish Balay   t[i++] = 1.7500;
219*9371c9d4SSatish Balay   y[i]   = 22.1250;
220*9371c9d4SSatish Balay   t[i++] = 2.2500;
221*9371c9d4SSatish Balay   y[i]   = 17.5125;
222*9371c9d4SSatish Balay   t[i++] = 2.7500;
223*9371c9d4SSatish Balay   y[i]   = 14.2500;
224*9371c9d4SSatish Balay   t[i++] = 3.2500;
225*9371c9d4SSatish Balay   y[i]   = 9.4500;
226*9371c9d4SSatish Balay   t[i++] = 3.7500;
227*9371c9d4SSatish Balay   y[i]   = 9.1500;
228*9371c9d4SSatish Balay   t[i++] = 4.2500;
229*9371c9d4SSatish Balay   y[i]   = 7.9125;
230*9371c9d4SSatish Balay   t[i++] = 4.7500;
231*9371c9d4SSatish Balay   y[i]   = 8.4750;
232*9371c9d4SSatish Balay   t[i++] = 5.2500;
233*9371c9d4SSatish Balay   y[i]   = 6.1125;
234*9371c9d4SSatish Balay   t[i++] = 5.7500;
235*9371c9d4SSatish Balay   y[i]   = 80.0000;
236*9371c9d4SSatish Balay   t[i++] = 0.5000;
237*9371c9d4SSatish Balay   y[i]   = 79.0000;
238*9371c9d4SSatish Balay   t[i++] = 0.6250;
239*9371c9d4SSatish Balay   y[i]   = 63.8000;
240*9371c9d4SSatish Balay   t[i++] = 0.7500;
241*9371c9d4SSatish Balay   y[i]   = 57.2000;
242*9371c9d4SSatish Balay   t[i++] = 0.8750;
243*9371c9d4SSatish Balay   y[i]   = 53.2000;
244*9371c9d4SSatish Balay   t[i++] = 1.0000;
245*9371c9d4SSatish Balay   y[i]   = 42.5000;
246*9371c9d4SSatish Balay   t[i++] = 1.2500;
247*9371c9d4SSatish Balay   y[i]   = 26.8000;
248*9371c9d4SSatish Balay   t[i++] = 1.7500;
249*9371c9d4SSatish Balay   y[i]   = 20.4000;
250*9371c9d4SSatish Balay   t[i++] = 2.2500;
251*9371c9d4SSatish Balay   y[i]   = 26.8500;
252*9371c9d4SSatish Balay   t[i++] = 1.7500;
253*9371c9d4SSatish Balay   y[i]   = 21.0000;
254*9371c9d4SSatish Balay   t[i++] = 2.2500;
255*9371c9d4SSatish Balay   y[i]   = 16.4625;
256*9371c9d4SSatish Balay   t[i++] = 2.7500;
257*9371c9d4SSatish Balay   y[i]   = 12.5250;
258*9371c9d4SSatish Balay   t[i++] = 3.2500;
259*9371c9d4SSatish Balay   y[i]   = 10.5375;
260*9371c9d4SSatish Balay   t[i++] = 3.7500;
261*9371c9d4SSatish Balay   y[i]   = 8.5875;
262*9371c9d4SSatish Balay   t[i++] = 4.2500;
263*9371c9d4SSatish Balay   y[i]   = 7.1250;
264*9371c9d4SSatish Balay   t[i++] = 4.7500;
265*9371c9d4SSatish Balay   y[i]   = 6.1125;
266*9371c9d4SSatish Balay   t[i++] = 5.2500;
267*9371c9d4SSatish Balay   y[i]   = 5.9625;
268*9371c9d4SSatish Balay   t[i++] = 5.7500;
269*9371c9d4SSatish Balay   y[i]   = 74.1000;
270*9371c9d4SSatish Balay   t[i++] = 0.5000;
271*9371c9d4SSatish Balay   y[i]   = 67.3000;
272*9371c9d4SSatish Balay   t[i++] = 0.6250;
273*9371c9d4SSatish Balay   y[i]   = 60.8000;
274*9371c9d4SSatish Balay   t[i++] = 0.7500;
275*9371c9d4SSatish Balay   y[i]   = 55.5000;
276*9371c9d4SSatish Balay   t[i++] = 0.8750;
277*9371c9d4SSatish Balay   y[i]   = 50.3000;
278*9371c9d4SSatish Balay   t[i++] = 1.0000;
279*9371c9d4SSatish Balay   y[i]   = 41.0000;
280*9371c9d4SSatish Balay   t[i++] = 1.2500;
281*9371c9d4SSatish Balay   y[i]   = 29.4000;
282*9371c9d4SSatish Balay   t[i++] = 1.7500;
283*9371c9d4SSatish Balay   y[i]   = 20.4000;
284*9371c9d4SSatish Balay   t[i++] = 2.2500;
285*9371c9d4SSatish Balay   y[i]   = 29.3625;
286*9371c9d4SSatish Balay   t[i++] = 1.7500;
287*9371c9d4SSatish Balay   y[i]   = 21.1500;
288*9371c9d4SSatish Balay   t[i++] = 2.2500;
289*9371c9d4SSatish Balay   y[i]   = 16.7625;
290*9371c9d4SSatish Balay   t[i++] = 2.7500;
291*9371c9d4SSatish Balay   y[i]   = 13.2000;
292*9371c9d4SSatish Balay   t[i++] = 3.2500;
293*9371c9d4SSatish Balay   y[i]   = 10.8750;
294*9371c9d4SSatish Balay   t[i++] = 3.7500;
295*9371c9d4SSatish Balay   y[i]   = 8.1750;
296*9371c9d4SSatish Balay   t[i++] = 4.2500;
297*9371c9d4SSatish Balay   y[i]   = 7.3500;
298*9371c9d4SSatish Balay   t[i++] = 4.7500;
299*9371c9d4SSatish Balay   y[i]   = 5.9625;
300*9371c9d4SSatish Balay   t[i++] = 5.2500;
301*9371c9d4SSatish Balay   y[i]   = 5.6250;
302*9371c9d4SSatish Balay   t[i++] = 5.7500;
303*9371c9d4SSatish Balay   y[i]   = 81.5000;
304*9371c9d4SSatish Balay   t[i++] = .5000;
305*9371c9d4SSatish Balay   y[i]   = 62.4000;
306*9371c9d4SSatish Balay   t[i++] = .7500;
307*9371c9d4SSatish Balay   y[i]   = 32.5000;
308*9371c9d4SSatish Balay   t[i++] = 1.5000;
309*9371c9d4SSatish Balay   y[i]   = 12.4100;
310*9371c9d4SSatish Balay   t[i++] = 3.0000;
311*9371c9d4SSatish Balay   y[i]   = 13.1200;
312*9371c9d4SSatish Balay   t[i++] = 3.0000;
313*9371c9d4SSatish Balay   y[i]   = 15.5600;
314*9371c9d4SSatish Balay   t[i++] = 3.0000;
315*9371c9d4SSatish Balay   y[i]   = 5.6300;
316*9371c9d4SSatish Balay   t[i++] = 6.0000;
317*9371c9d4SSatish Balay   y[i]   = 78.0000;
318*9371c9d4SSatish Balay   t[i++] = .5000;
319*9371c9d4SSatish Balay   y[i]   = 59.9000;
320*9371c9d4SSatish Balay   t[i++] = .7500;
321*9371c9d4SSatish Balay   y[i]   = 33.2000;
322*9371c9d4SSatish Balay   t[i++] = 1.5000;
323*9371c9d4SSatish Balay   y[i]   = 13.8400;
324*9371c9d4SSatish Balay   t[i++] = 3.0000;
325*9371c9d4SSatish Balay   y[i]   = 12.7500;
326*9371c9d4SSatish Balay   t[i++] = 3.0000;
327*9371c9d4SSatish Balay   y[i]   = 14.6200;
328*9371c9d4SSatish Balay   t[i++] = 3.0000;
329*9371c9d4SSatish Balay   y[i]   = 3.9400;
330*9371c9d4SSatish Balay   t[i++] = 6.0000;
331*9371c9d4SSatish Balay   y[i]   = 76.8000;
332*9371c9d4SSatish Balay   t[i++] = .5000;
333*9371c9d4SSatish Balay   y[i]   = 61.0000;
334*9371c9d4SSatish Balay   t[i++] = .7500;
335*9371c9d4SSatish Balay   y[i]   = 32.9000;
336*9371c9d4SSatish Balay   t[i++] = 1.5000;
337*9371c9d4SSatish Balay   y[i]   = 13.8700;
338*9371c9d4SSatish Balay   t[i++] = 3.0000;
339*9371c9d4SSatish Balay   y[i]   = 11.8100;
340*9371c9d4SSatish Balay   t[i++] = 3.0000;
341*9371c9d4SSatish Balay   y[i]   = 13.3100;
342*9371c9d4SSatish Balay   t[i++] = 3.0000;
343*9371c9d4SSatish Balay   y[i]   = 5.4400;
344*9371c9d4SSatish Balay   t[i++] = 6.0000;
345*9371c9d4SSatish Balay   y[i]   = 78.0000;
346*9371c9d4SSatish Balay   t[i++] = .5000;
347*9371c9d4SSatish Balay   y[i]   = 63.5000;
348*9371c9d4SSatish Balay   t[i++] = .7500;
349*9371c9d4SSatish Balay   y[i]   = 33.8000;
350*9371c9d4SSatish Balay   t[i++] = 1.5000;
351*9371c9d4SSatish Balay   y[i]   = 12.5600;
352*9371c9d4SSatish Balay   t[i++] = 3.0000;
353*9371c9d4SSatish Balay   y[i]   = 5.6300;
354*9371c9d4SSatish Balay   t[i++] = 6.0000;
355*9371c9d4SSatish Balay   y[i]   = 12.7500;
356*9371c9d4SSatish Balay   t[i++] = 3.0000;
357*9371c9d4SSatish Balay   y[i]   = 13.1200;
358*9371c9d4SSatish Balay   t[i++] = 3.0000;
359*9371c9d4SSatish Balay   y[i]   = 5.4400;
360*9371c9d4SSatish Balay   t[i++] = 6.0000;
361*9371c9d4SSatish Balay   y[i]   = 76.8000;
362*9371c9d4SSatish Balay   t[i++] = .5000;
363*9371c9d4SSatish Balay   y[i]   = 60.0000;
364*9371c9d4SSatish Balay   t[i++] = .7500;
365*9371c9d4SSatish Balay   y[i]   = 47.8000;
366*9371c9d4SSatish Balay   t[i++] = 1.0000;
367*9371c9d4SSatish Balay   y[i]   = 32.0000;
368*9371c9d4SSatish Balay   t[i++] = 1.5000;
369*9371c9d4SSatish Balay   y[i]   = 22.2000;
370*9371c9d4SSatish Balay   t[i++] = 2.0000;
371*9371c9d4SSatish Balay   y[i]   = 22.5700;
372*9371c9d4SSatish Balay   t[i++] = 2.0000;
373*9371c9d4SSatish Balay   y[i]   = 18.8200;
374*9371c9d4SSatish Balay   t[i++] = 2.5000;
375*9371c9d4SSatish Balay   y[i]   = 13.9500;
376*9371c9d4SSatish Balay   t[i++] = 3.0000;
377*9371c9d4SSatish Balay   y[i]   = 11.2500;
378*9371c9d4SSatish Balay   t[i++] = 4.0000;
379*9371c9d4SSatish Balay   y[i]   = 9.0000;
380*9371c9d4SSatish Balay   t[i++] = 5.0000;
381*9371c9d4SSatish Balay   y[i]   = 6.6700;
382*9371c9d4SSatish Balay   t[i++] = 6.0000;
383*9371c9d4SSatish Balay   y[i]   = 75.8000;
384*9371c9d4SSatish Balay   t[i++] = .5000;
385*9371c9d4SSatish Balay   y[i]   = 62.0000;
386*9371c9d4SSatish Balay   t[i++] = .7500;
387*9371c9d4SSatish Balay   y[i]   = 48.8000;
388*9371c9d4SSatish Balay   t[i++] = 1.0000;
389*9371c9d4SSatish Balay   y[i]   = 35.2000;
390*9371c9d4SSatish Balay   t[i++] = 1.5000;
391*9371c9d4SSatish Balay   y[i]   = 20.0000;
392*9371c9d4SSatish Balay   t[i++] = 2.0000;
393*9371c9d4SSatish Balay   y[i]   = 20.3200;
394*9371c9d4SSatish Balay   t[i++] = 2.0000;
395*9371c9d4SSatish Balay   y[i]   = 19.3100;
396*9371c9d4SSatish Balay   t[i++] = 2.5000;
397*9371c9d4SSatish Balay   y[i]   = 12.7500;
398*9371c9d4SSatish Balay   t[i++] = 3.0000;
399*9371c9d4SSatish Balay   y[i]   = 10.4200;
400*9371c9d4SSatish Balay   t[i++] = 4.0000;
401*9371c9d4SSatish Balay   y[i]   = 7.3100;
402*9371c9d4SSatish Balay   t[i++] = 5.0000;
403*9371c9d4SSatish Balay   y[i]   = 7.4200;
404*9371c9d4SSatish Balay   t[i++] = 6.0000;
405*9371c9d4SSatish Balay   y[i]   = 70.5000;
406*9371c9d4SSatish Balay   t[i++] = .5000;
407*9371c9d4SSatish Balay   y[i]   = 59.5000;
408*9371c9d4SSatish Balay   t[i++] = .7500;
409*9371c9d4SSatish Balay   y[i]   = 48.5000;
410*9371c9d4SSatish Balay   t[i++] = 1.0000;
411*9371c9d4SSatish Balay   y[i]   = 35.8000;
412*9371c9d4SSatish Balay   t[i++] = 1.5000;
413*9371c9d4SSatish Balay   y[i]   = 21.0000;
414*9371c9d4SSatish Balay   t[i++] = 2.0000;
415*9371c9d4SSatish Balay   y[i]   = 21.6700;
416*9371c9d4SSatish Balay   t[i++] = 2.0000;
417*9371c9d4SSatish Balay   y[i]   = 21.0000;
418*9371c9d4SSatish Balay   t[i++] = 2.5000;
419*9371c9d4SSatish Balay   y[i]   = 15.6400;
420*9371c9d4SSatish Balay   t[i++] = 3.0000;
421*9371c9d4SSatish Balay   y[i]   = 8.1700;
422*9371c9d4SSatish Balay   t[i++] = 4.0000;
423*9371c9d4SSatish Balay   y[i]   = 8.5500;
424*9371c9d4SSatish Balay   t[i++] = 5.0000;
425*9371c9d4SSatish Balay   y[i]   = 10.1200;
426*9371c9d4SSatish Balay   t[i++] = 6.0000;
427*9371c9d4SSatish Balay   y[i]   = 78.0000;
428*9371c9d4SSatish Balay   t[i++] = .5000;
429*9371c9d4SSatish Balay   y[i]   = 66.0000;
430*9371c9d4SSatish Balay   t[i++] = .6250;
431*9371c9d4SSatish Balay   y[i]   = 62.0000;
432*9371c9d4SSatish Balay   t[i++] = .7500;
433*9371c9d4SSatish Balay   y[i]   = 58.0000;
434*9371c9d4SSatish Balay   t[i++] = .8750;
435*9371c9d4SSatish Balay   y[i]   = 47.7000;
436*9371c9d4SSatish Balay   t[i++] = 1.0000;
437*9371c9d4SSatish Balay   y[i]   = 37.8000;
438*9371c9d4SSatish Balay   t[i++] = 1.2500;
439*9371c9d4SSatish Balay   y[i]   = 20.2000;
440*9371c9d4SSatish Balay   t[i++] = 2.2500;
441*9371c9d4SSatish Balay   y[i]   = 21.0700;
442*9371c9d4SSatish Balay   t[i++] = 2.2500;
443*9371c9d4SSatish Balay   y[i]   = 13.8700;
444*9371c9d4SSatish Balay   t[i++] = 2.7500;
445*9371c9d4SSatish Balay   y[i]   = 9.6700;
446*9371c9d4SSatish Balay   t[i++] = 3.2500;
447*9371c9d4SSatish Balay   y[i]   = 7.7600;
448*9371c9d4SSatish Balay   t[i++] = 3.7500;
449*9371c9d4SSatish Balay   y[i]   = 5.4400;
450*9371c9d4SSatish Balay   t[i++] = 4.2500;
451*9371c9d4SSatish Balay   y[i]   = 4.8700;
452*9371c9d4SSatish Balay   t[i++] = 4.7500;
453*9371c9d4SSatish Balay   y[i]   = 4.0100;
454*9371c9d4SSatish Balay   t[i++] = 5.2500;
455*9371c9d4SSatish Balay   y[i]   = 3.7500;
456*9371c9d4SSatish Balay   t[i++] = 5.7500;
457*9371c9d4SSatish Balay   y[i]   = 24.1900;
458*9371c9d4SSatish Balay   t[i++] = 3.0000;
459*9371c9d4SSatish Balay   y[i]   = 25.7600;
460*9371c9d4SSatish Balay   t[i++] = 3.0000;
461*9371c9d4SSatish Balay   y[i]   = 18.0700;
462*9371c9d4SSatish Balay   t[i++] = 3.0000;
463*9371c9d4SSatish Balay   y[i]   = 11.8100;
464*9371c9d4SSatish Balay   t[i++] = 3.0000;
465*9371c9d4SSatish Balay   y[i]   = 12.0700;
466*9371c9d4SSatish Balay   t[i++] = 3.0000;
467*9371c9d4SSatish Balay   y[i]   = 16.1200;
468*9371c9d4SSatish Balay   t[i++] = 3.0000;
469*9371c9d4SSatish Balay   y[i]   = 70.8000;
470*9371c9d4SSatish Balay   t[i++] = .5000;
471*9371c9d4SSatish Balay   y[i]   = 54.7000;
472*9371c9d4SSatish Balay   t[i++] = .7500;
473*9371c9d4SSatish Balay   y[i]   = 48.0000;
474*9371c9d4SSatish Balay   t[i++] = 1.0000;
475*9371c9d4SSatish Balay   y[i]   = 39.8000;
476*9371c9d4SSatish Balay   t[i++] = 1.5000;
477*9371c9d4SSatish Balay   y[i]   = 29.8000;
478*9371c9d4SSatish Balay   t[i++] = 2.0000;
479*9371c9d4SSatish Balay   y[i]   = 23.7000;
480*9371c9d4SSatish Balay   t[i++] = 2.5000;
481*9371c9d4SSatish Balay   y[i]   = 29.6200;
482*9371c9d4SSatish Balay   t[i++] = 2.0000;
483*9371c9d4SSatish Balay   y[i]   = 23.8100;
484*9371c9d4SSatish Balay   t[i++] = 2.5000;
485*9371c9d4SSatish Balay   y[i]   = 17.7000;
486*9371c9d4SSatish Balay   t[i++] = 3.0000;
487*9371c9d4SSatish Balay   y[i]   = 11.5500;
488*9371c9d4SSatish Balay   t[i++] = 4.0000;
489*9371c9d4SSatish Balay   y[i]   = 12.0700;
490*9371c9d4SSatish Balay   t[i++] = 5.0000;
491*9371c9d4SSatish Balay   y[i]   = 8.7400;
492*9371c9d4SSatish Balay   t[i++] = 6.0000;
493*9371c9d4SSatish Balay   y[i]   = 80.7000;
494*9371c9d4SSatish Balay   t[i++] = .5000;
495*9371c9d4SSatish Balay   y[i]   = 61.3000;
496*9371c9d4SSatish Balay   t[i++] = .7500;
497*9371c9d4SSatish Balay   y[i]   = 47.5000;
498*9371c9d4SSatish Balay   t[i++] = 1.0000;
499*9371c9d4SSatish Balay   y[i]   = 29.0000;
500*9371c9d4SSatish Balay   t[i++] = 1.5000;
501*9371c9d4SSatish Balay   y[i]   = 24.0000;
502*9371c9d4SSatish Balay   t[i++] = 2.0000;
503*9371c9d4SSatish Balay   y[i]   = 17.7000;
504*9371c9d4SSatish Balay   t[i++] = 2.5000;
505*9371c9d4SSatish Balay   y[i]   = 24.5600;
506*9371c9d4SSatish Balay   t[i++] = 2.0000;
507*9371c9d4SSatish Balay   y[i]   = 18.6700;
508*9371c9d4SSatish Balay   t[i++] = 2.5000;
509*9371c9d4SSatish Balay   y[i]   = 16.2400;
510*9371c9d4SSatish Balay   t[i++] = 3.0000;
511*9371c9d4SSatish Balay   y[i]   = 8.7400;
512*9371c9d4SSatish Balay   t[i++] = 4.0000;
513*9371c9d4SSatish Balay   y[i]   = 7.8700;
514*9371c9d4SSatish Balay   t[i++] = 5.0000;
515*9371c9d4SSatish Balay   y[i]   = 8.5100;
516*9371c9d4SSatish Balay   t[i++] = 6.0000;
517*9371c9d4SSatish Balay   y[i]   = 66.7000;
518*9371c9d4SSatish Balay   t[i++] = .5000;
519*9371c9d4SSatish Balay   y[i]   = 59.2000;
520*9371c9d4SSatish Balay   t[i++] = .7500;
521*9371c9d4SSatish Balay   y[i]   = 40.8000;
522*9371c9d4SSatish Balay   t[i++] = 1.0000;
523*9371c9d4SSatish Balay   y[i]   = 30.7000;
524*9371c9d4SSatish Balay   t[i++] = 1.5000;
525*9371c9d4SSatish Balay   y[i]   = 25.7000;
526*9371c9d4SSatish Balay   t[i++] = 2.0000;
527*9371c9d4SSatish Balay   y[i]   = 16.3000;
528*9371c9d4SSatish Balay   t[i++] = 2.5000;
529*9371c9d4SSatish Balay   y[i]   = 25.9900;
530*9371c9d4SSatish Balay   t[i++] = 2.0000;
531*9371c9d4SSatish Balay   y[i]   = 16.9500;
532*9371c9d4SSatish Balay   t[i++] = 2.5000;
533*9371c9d4SSatish Balay   y[i]   = 13.3500;
534*9371c9d4SSatish Balay   t[i++] = 3.0000;
535*9371c9d4SSatish Balay   y[i]   = 8.6200;
536*9371c9d4SSatish Balay   t[i++] = 4.0000;
537*9371c9d4SSatish Balay   y[i]   = 7.2000;
538*9371c9d4SSatish Balay   t[i++] = 5.0000;
539*9371c9d4SSatish Balay   y[i]   = 6.6400;
540*9371c9d4SSatish Balay   t[i++] = 6.0000;
541*9371c9d4SSatish Balay   y[i]   = 13.6900;
542*9371c9d4SSatish Balay   t[i++] = 3.0000;
543*9371c9d4SSatish Balay   y[i]   = 81.0000;
544*9371c9d4SSatish Balay   t[i++] = .5000;
545*9371c9d4SSatish Balay   y[i]   = 64.5000;
546*9371c9d4SSatish Balay   t[i++] = .7500;
547*9371c9d4SSatish Balay   y[i]   = 35.5000;
548*9371c9d4SSatish Balay   t[i++] = 1.5000;
549*9371c9d4SSatish Balay   y[i]   = 13.3100;
550*9371c9d4SSatish Balay   t[i++] = 3.0000;
551*9371c9d4SSatish Balay   y[i]   = 4.8700;
552*9371c9d4SSatish Balay   t[i++] = 6.0000;
553*9371c9d4SSatish Balay   y[i]   = 12.9400;
554*9371c9d4SSatish Balay   t[i++] = 3.0000;
555*9371c9d4SSatish Balay   y[i]   = 5.0600;
556*9371c9d4SSatish Balay   t[i++] = 6.0000;
557*9371c9d4SSatish Balay   y[i]   = 15.1900;
558*9371c9d4SSatish Balay   t[i++] = 3.0000;
559*9371c9d4SSatish Balay   y[i]   = 14.6200;
560*9371c9d4SSatish Balay   t[i++] = 3.0000;
561*9371c9d4SSatish Balay   y[i]   = 15.6400;
562*9371c9d4SSatish Balay   t[i++] = 3.0000;
563*9371c9d4SSatish Balay   y[i]   = 25.5000;
564*9371c9d4SSatish Balay   t[i++] = 1.7500;
565*9371c9d4SSatish Balay   y[i]   = 25.9500;
566*9371c9d4SSatish Balay   t[i++] = 1.7500;
567*9371c9d4SSatish Balay   y[i]   = 81.7000;
568*9371c9d4SSatish Balay   t[i++] = .5000;
569*9371c9d4SSatish Balay   y[i]   = 61.6000;
570*9371c9d4SSatish Balay   t[i++] = .7500;
571*9371c9d4SSatish Balay   y[i]   = 29.8000;
572*9371c9d4SSatish Balay   t[i++] = 1.7500;
573*9371c9d4SSatish Balay   y[i]   = 29.8100;
574*9371c9d4SSatish Balay   t[i++] = 1.7500;
575*9371c9d4SSatish Balay   y[i]   = 17.1700;
576*9371c9d4SSatish Balay   t[i++] = 2.7500;
577*9371c9d4SSatish Balay   y[i]   = 10.3900;
578*9371c9d4SSatish Balay   t[i++] = 3.7500;
579*9371c9d4SSatish Balay   y[i]   = 28.4000;
580*9371c9d4SSatish Balay   t[i++] = 1.7500;
581*9371c9d4SSatish Balay   y[i]   = 28.6900;
582*9371c9d4SSatish Balay   t[i++] = 1.7500;
583*9371c9d4SSatish Balay   y[i]   = 81.3000;
584*9371c9d4SSatish Balay   t[i++] = .5000;
585*9371c9d4SSatish Balay   y[i]   = 60.9000;
586*9371c9d4SSatish Balay   t[i++] = .7500;
587*9371c9d4SSatish Balay   y[i]   = 16.6500;
588*9371c9d4SSatish Balay   t[i++] = 2.7500;
589*9371c9d4SSatish Balay   y[i]   = 10.0500;
590*9371c9d4SSatish Balay   t[i++] = 3.7500;
591*9371c9d4SSatish Balay   y[i]   = 28.9000;
592*9371c9d4SSatish Balay   t[i++] = 1.7500;
593*9371c9d4SSatish Balay   y[i]   = 28.9500;
594*9371c9d4SSatish Balay   t[i++] = 1.7500;
595c4762a1bSJed Brown   PetscFunctionReturn(0);
596c4762a1bSJed Brown }
597c4762a1bSJed Brown 
598c4762a1bSJed Brown /*TEST
599c4762a1bSJed Brown 
600c4762a1bSJed Brown    build:
601c4762a1bSJed Brown       requires: !complex !single
602c4762a1bSJed Brown 
603c4762a1bSJed Brown    test:
604c4762a1bSJed Brown       args: -tao_smonitor -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-5
605c4762a1bSJed Brown 
606c4762a1bSJed Brown    test:
607c4762a1bSJed Brown       suffix: 2
608c4762a1bSJed 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
609c4762a1bSJed Brown 
610c4762a1bSJed Brown    test:
611c4762a1bSJed Brown       suffix: 3
612c4762a1bSJed 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
613c4762a1bSJed Brown 
614cd1c4666STristan Konolige    test:
615cd1c4666STristan Konolige       suffix: 4
616cd1c4666STristan 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
617cd1c4666STristan Konolige 
618c4762a1bSJed Brown TEST*/
619