xref: /petsc/src/tao/leastsquares/tests/chwirut2.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
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  This version tests correlated terms using both vector and listed forms
10c4762a1bSJed Brown */
11c4762a1bSJed Brown 
12c4762a1bSJed Brown #include <petsctao.h>
13c4762a1bSJed Brown 
14c4762a1bSJed Brown /*
15c4762a1bSJed Brown Description:   These data are the result of a NIST study involving
16c4762a1bSJed Brown                ultrasonic calibration.  The response variable is
17c4762a1bSJed Brown                ultrasonic response, and the predictor variable is
18c4762a1bSJed Brown                metal distance.
19c4762a1bSJed Brown 
20c4762a1bSJed Brown Reference:     Chwirut, D., NIST (197?).
21c4762a1bSJed Brown                Ultrasonic Reference Block Study.
22c4762a1bSJed Brown */
23c4762a1bSJed Brown 
24c4762a1bSJed Brown static char help[] = "Finds the nonlinear least-squares solution to the model \n\
25c4762a1bSJed Brown             y = exp[-b1*x]/(b2+b3*x)  +  e \n";
26c4762a1bSJed Brown 
27c4762a1bSJed Brown #define NOBSERVATIONS 214
28c4762a1bSJed Brown #define NPARAMETERS   3
29c4762a1bSJed Brown 
30c4762a1bSJed Brown /* User-defined application context */
31c4762a1bSJed Brown typedef struct {
32c4762a1bSJed Brown   /* Working space */
33c4762a1bSJed Brown   PetscReal t[NOBSERVATIONS];              /* array of independent variables of observation */
34c4762a1bSJed Brown   PetscReal y[NOBSERVATIONS];              /* array of dependent variables */
35c4762a1bSJed Brown   PetscReal j[NOBSERVATIONS][NPARAMETERS]; /* dense jacobian matrix array*/
36c4762a1bSJed Brown   PetscInt  idm[NOBSERVATIONS];            /* Matrix indices for jacobian */
37c4762a1bSJed Brown   PetscInt  idn[NPARAMETERS];
38c4762a1bSJed Brown } AppCtx;
39c4762a1bSJed Brown 
40c4762a1bSJed Brown /* User provided Routines */
41c4762a1bSJed Brown PetscErrorCode InitializeData(AppCtx *user);
42c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec);
43c4762a1bSJed Brown PetscErrorCode EvaluateFunction(Tao, Vec, Vec, void *);
44c4762a1bSJed Brown PetscErrorCode EvaluateJacobian(Tao, Vec, Mat, Mat, void *);
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /*--------------------------------------------------------------------*/
main(int argc,char ** argv)47d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
48d71ae5a4SJacob Faibussowitsch {
49c4762a1bSJed Brown   PetscInt  wtype = 0;
50c4762a1bSJed Brown   Vec       x, f; /* solution, function */
51c4762a1bSJed Brown   Vec       w;    /* weights */
52c4762a1bSJed Brown   Mat       J;    /* Jacobian matrix */
53c4762a1bSJed Brown   Tao       tao;  /* Tao solver context */
54c4762a1bSJed Brown   PetscInt  i;    /* iteration information */
55c4762a1bSJed Brown   PetscReal hist[100], resid[100];
56c4762a1bSJed Brown   PetscInt  lits[100];
57c4762a1bSJed Brown   PetscInt  w_row[NOBSERVATIONS]; /* explicit weights */
58c4762a1bSJed Brown   PetscInt  w_col[NOBSERVATIONS];
59c4762a1bSJed Brown   PetscReal w_vals[NOBSERVATIONS];
60c4762a1bSJed Brown   PetscBool flg;
61c4762a1bSJed Brown   AppCtx    user; /* user-defined work context */
62c4762a1bSJed Brown 
63327415f7SBarry Smith   PetscFunctionBeginUser;
64*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-wtype", &wtype, &flg));
6663a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "wtype=%" PetscInt_FMT "\n", wtype));
67c4762a1bSJed Brown   /* Allocate vectors */
689566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(MPI_COMM_SELF, NPARAMETERS, &x));
699566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(MPI_COMM_SELF, NOBSERVATIONS, &f));
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(f, &w));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* no correlation, but set in different ways */
749566063dSJacob Faibussowitsch   PetscCall(VecSet(w, 1.0));
75c4762a1bSJed Brown   for (i = 0; i < NOBSERVATIONS; i++) {
769371c9d4SSatish Balay     w_row[i]  = i;
779371c9d4SSatish Balay     w_col[i]  = i;
789371c9d4SSatish Balay     w_vals[i] = 1.0;
79c4762a1bSJed Brown   }
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* Create the Jacobian matrix. */
829566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(MPI_COMM_SELF, NOBSERVATIONS, NPARAMETERS, NULL, &J));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   for (i = 0; i < NOBSERVATIONS; i++) user.idm[i] = i;
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   for (i = 0; i < NPARAMETERS; i++) user.idn[i] = i;
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
899566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_SELF, &tao));
909566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOPOUNDERS));
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   /* Set the function and Jacobian routines. */
939566063dSJacob Faibussowitsch   PetscCall(InitializeData(&user));
949566063dSJacob Faibussowitsch   PetscCall(FormStartingPoint(x));
959566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, x));
969566063dSJacob Faibussowitsch   PetscCall(TaoSetResidualRoutine(tao, f, EvaluateFunction, (void *)&user));
97c4762a1bSJed Brown   if (wtype == 1) {
989566063dSJacob Faibussowitsch     PetscCall(TaoSetResidualWeights(tao, w, 0, NULL, NULL, NULL));
99c4762a1bSJed Brown   } else if (wtype == 2) {
1009566063dSJacob Faibussowitsch     PetscCall(TaoSetResidualWeights(tao, NULL, NOBSERVATIONS, w_row, w_col, w_vals));
101c4762a1bSJed Brown   }
1029566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianResidualRoutine(tao, J, J, EvaluateJacobian, (void *)&user));
103606f75f6SBarry Smith   PetscCall(TaoSetTolerances(tao, 1e-5, 0.0, PETSC_CURRENT));
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* Check for any TAO command line arguments */
1069566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
107c4762a1bSJed Brown 
1089566063dSJacob Faibussowitsch   PetscCall(TaoSetConvergenceHistory(tao, hist, resid, 0, lits, 100, PETSC_TRUE));
109c4762a1bSJed Brown   /* Perform the Solve */
1109566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /* Free TAO data structures */
1139566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* Free PETSc data structures */
1169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&w));
1189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&f));
1199566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
122b122ec5aSJacob Faibussowitsch   return 0;
123c4762a1bSJed Brown }
124c4762a1bSJed Brown 
125c4762a1bSJed Brown /*--------------------------------------------------------------------*/
EvaluateFunction(Tao tao,Vec X,Vec F,void * ptr)126d71ae5a4SJacob Faibussowitsch PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr)
127d71ae5a4SJacob Faibussowitsch {
128c4762a1bSJed Brown   AppCtx          *user = (AppCtx *)ptr;
129c4762a1bSJed Brown   PetscInt         i;
13005579b36STristan Konolige   PetscReal       *y = user->y, *f, *t = user->t;
13105579b36STristan Konolige   const PetscReal *x;
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   PetscFunctionBegin;
1349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1359566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
136c4762a1bSJed Brown 
137ad540459SPierre Jolivet   for (i = 0; i < NOBSERVATIONS; i++) f[i] = y[i] - PetscExpScalar(-x[0] * t[i]) / (x[1] + x[2] * t[i]);
1389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1403ba16761SJacob Faibussowitsch   PetscCall(PetscLogFlops(6 * NOBSERVATIONS));
1413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
142c4762a1bSJed Brown }
143c4762a1bSJed Brown 
144c4762a1bSJed Brown /*------------------------------------------------------------*/
145c4762a1bSJed Brown /* J[i][j] = df[i]/dt[j] */
EvaluateJacobian(Tao tao,Vec X,Mat J,Mat Jpre,void * ptr)146d71ae5a4SJacob Faibussowitsch PetscErrorCode EvaluateJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr)
147d71ae5a4SJacob Faibussowitsch {
148c4762a1bSJed Brown   AppCtx          *user = (AppCtx *)ptr;
149c4762a1bSJed Brown   PetscInt         i;
15005579b36STristan Konolige   PetscReal       *t = user->t;
15105579b36STristan Konolige   const PetscReal *x;
152c4762a1bSJed Brown   PetscReal        base;
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   PetscFunctionBegin;
1559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
156c4762a1bSJed Brown   for (i = 0; i < NOBSERVATIONS; i++) {
157c4762a1bSJed Brown     base = PetscExpScalar(-x[0] * t[i]) / (x[1] + x[2] * t[i]);
158c4762a1bSJed Brown 
159c4762a1bSJed Brown     user->j[i][0] = t[i] * base;
160c4762a1bSJed Brown     user->j[i][1] = base / (x[1] + x[2] * t[i]);
161c4762a1bSJed Brown     user->j[i][2] = base * t[i] / (x[1] + x[2] * t[i]);
162c4762a1bSJed Brown   }
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* Assemble the matrix */
1659566063dSJacob Faibussowitsch   PetscCall(MatSetValues(J, NOBSERVATIONS, user->idm, NPARAMETERS, user->idn, (PetscReal *)user->j, INSERT_VALUES));
1669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
168c4762a1bSJed Brown 
1699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1703ba16761SJacob Faibussowitsch   PetscCall(PetscLogFlops(NOBSERVATIONS * 13));
1713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
172c4762a1bSJed Brown }
173c4762a1bSJed Brown 
174c4762a1bSJed Brown /* ------------------------------------------------------------ */
FormStartingPoint(Vec X)175d71ae5a4SJacob Faibussowitsch PetscErrorCode FormStartingPoint(Vec X)
176d71ae5a4SJacob Faibussowitsch {
177c4762a1bSJed Brown   PetscReal *x;
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   PetscFunctionBegin;
1809566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
181c4762a1bSJed Brown   x[0] = 1.19;
182c4762a1bSJed Brown   x[1] = -1.86;
183c4762a1bSJed Brown   x[2] = 1.08;
1849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
1853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
186c4762a1bSJed Brown }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown /* ---------------------------------------------------------------------- */
InitializeData(AppCtx * user)189d71ae5a4SJacob Faibussowitsch PetscErrorCode InitializeData(AppCtx *user)
190d71ae5a4SJacob Faibussowitsch {
191c4762a1bSJed Brown   PetscReal *t = user->t, *y = user->y;
192c4762a1bSJed Brown   PetscInt   i = 0;
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   PetscFunctionBegin;
1959371c9d4SSatish Balay   y[i]   = 92.9000;
1969371c9d4SSatish Balay   t[i++] = 0.5000;
1979371c9d4SSatish Balay   y[i]   = 78.7000;
1989371c9d4SSatish Balay   t[i++] = 0.6250;
1999371c9d4SSatish Balay   y[i]   = 64.2000;
2009371c9d4SSatish Balay   t[i++] = 0.7500;
2019371c9d4SSatish Balay   y[i]   = 64.9000;
2029371c9d4SSatish Balay   t[i++] = 0.8750;
2039371c9d4SSatish Balay   y[i]   = 57.1000;
2049371c9d4SSatish Balay   t[i++] = 1.0000;
2059371c9d4SSatish Balay   y[i]   = 43.3000;
2069371c9d4SSatish Balay   t[i++] = 1.2500;
2079371c9d4SSatish Balay   y[i]   = 31.1000;
2089371c9d4SSatish Balay   t[i++] = 1.7500;
2099371c9d4SSatish Balay   y[i]   = 23.6000;
2109371c9d4SSatish Balay   t[i++] = 2.2500;
2119371c9d4SSatish Balay   y[i]   = 31.0500;
2129371c9d4SSatish Balay   t[i++] = 1.7500;
2139371c9d4SSatish Balay   y[i]   = 23.7750;
2149371c9d4SSatish Balay   t[i++] = 2.2500;
2159371c9d4SSatish Balay   y[i]   = 17.7375;
2169371c9d4SSatish Balay   t[i++] = 2.7500;
2179371c9d4SSatish Balay   y[i]   = 13.8000;
2189371c9d4SSatish Balay   t[i++] = 3.2500;
2199371c9d4SSatish Balay   y[i]   = 11.5875;
2209371c9d4SSatish Balay   t[i++] = 3.7500;
2219371c9d4SSatish Balay   y[i]   = 9.4125;
2229371c9d4SSatish Balay   t[i++] = 4.2500;
2239371c9d4SSatish Balay   y[i]   = 7.7250;
2249371c9d4SSatish Balay   t[i++] = 4.7500;
2259371c9d4SSatish Balay   y[i]   = 7.3500;
2269371c9d4SSatish Balay   t[i++] = 5.2500;
2279371c9d4SSatish Balay   y[i]   = 8.0250;
2289371c9d4SSatish Balay   t[i++] = 5.7500;
2299371c9d4SSatish Balay   y[i]   = 90.6000;
2309371c9d4SSatish Balay   t[i++] = 0.5000;
2319371c9d4SSatish Balay   y[i]   = 76.9000;
2329371c9d4SSatish Balay   t[i++] = 0.6250;
2339371c9d4SSatish Balay   y[i]   = 71.6000;
2349371c9d4SSatish Balay   t[i++] = 0.7500;
2359371c9d4SSatish Balay   y[i]   = 63.6000;
2369371c9d4SSatish Balay   t[i++] = 0.8750;
2379371c9d4SSatish Balay   y[i]   = 54.0000;
2389371c9d4SSatish Balay   t[i++] = 1.0000;
2399371c9d4SSatish Balay   y[i]   = 39.2000;
2409371c9d4SSatish Balay   t[i++] = 1.2500;
2419371c9d4SSatish Balay   y[i]   = 29.3000;
2429371c9d4SSatish Balay   t[i++] = 1.7500;
2439371c9d4SSatish Balay   y[i]   = 21.4000;
2449371c9d4SSatish Balay   t[i++] = 2.2500;
2459371c9d4SSatish Balay   y[i]   = 29.1750;
2469371c9d4SSatish Balay   t[i++] = 1.7500;
2479371c9d4SSatish Balay   y[i]   = 22.1250;
2489371c9d4SSatish Balay   t[i++] = 2.2500;
2499371c9d4SSatish Balay   y[i]   = 17.5125;
2509371c9d4SSatish Balay   t[i++] = 2.7500;
2519371c9d4SSatish Balay   y[i]   = 14.2500;
2529371c9d4SSatish Balay   t[i++] = 3.2500;
2539371c9d4SSatish Balay   y[i]   = 9.4500;
2549371c9d4SSatish Balay   t[i++] = 3.7500;
2559371c9d4SSatish Balay   y[i]   = 9.1500;
2569371c9d4SSatish Balay   t[i++] = 4.2500;
2579371c9d4SSatish Balay   y[i]   = 7.9125;
2589371c9d4SSatish Balay   t[i++] = 4.7500;
2599371c9d4SSatish Balay   y[i]   = 8.4750;
2609371c9d4SSatish Balay   t[i++] = 5.2500;
2619371c9d4SSatish Balay   y[i]   = 6.1125;
2629371c9d4SSatish Balay   t[i++] = 5.7500;
2639371c9d4SSatish Balay   y[i]   = 80.0000;
2649371c9d4SSatish Balay   t[i++] = 0.5000;
2659371c9d4SSatish Balay   y[i]   = 79.0000;
2669371c9d4SSatish Balay   t[i++] = 0.6250;
2679371c9d4SSatish Balay   y[i]   = 63.8000;
2689371c9d4SSatish Balay   t[i++] = 0.7500;
2699371c9d4SSatish Balay   y[i]   = 57.2000;
2709371c9d4SSatish Balay   t[i++] = 0.8750;
2719371c9d4SSatish Balay   y[i]   = 53.2000;
2729371c9d4SSatish Balay   t[i++] = 1.0000;
2739371c9d4SSatish Balay   y[i]   = 42.5000;
2749371c9d4SSatish Balay   t[i++] = 1.2500;
2759371c9d4SSatish Balay   y[i]   = 26.8000;
2769371c9d4SSatish Balay   t[i++] = 1.7500;
2779371c9d4SSatish Balay   y[i]   = 20.4000;
2789371c9d4SSatish Balay   t[i++] = 2.2500;
2799371c9d4SSatish Balay   y[i]   = 26.8500;
2809371c9d4SSatish Balay   t[i++] = 1.7500;
2819371c9d4SSatish Balay   y[i]   = 21.0000;
2829371c9d4SSatish Balay   t[i++] = 2.2500;
2839371c9d4SSatish Balay   y[i]   = 16.4625;
2849371c9d4SSatish Balay   t[i++] = 2.7500;
2859371c9d4SSatish Balay   y[i]   = 12.5250;
2869371c9d4SSatish Balay   t[i++] = 3.2500;
2879371c9d4SSatish Balay   y[i]   = 10.5375;
2889371c9d4SSatish Balay   t[i++] = 3.7500;
2899371c9d4SSatish Balay   y[i]   = 8.5875;
2909371c9d4SSatish Balay   t[i++] = 4.2500;
2919371c9d4SSatish Balay   y[i]   = 7.1250;
2929371c9d4SSatish Balay   t[i++] = 4.7500;
2939371c9d4SSatish Balay   y[i]   = 6.1125;
2949371c9d4SSatish Balay   t[i++] = 5.2500;
2959371c9d4SSatish Balay   y[i]   = 5.9625;
2969371c9d4SSatish Balay   t[i++] = 5.7500;
2979371c9d4SSatish Balay   y[i]   = 74.1000;
2989371c9d4SSatish Balay   t[i++] = 0.5000;
2999371c9d4SSatish Balay   y[i]   = 67.3000;
3009371c9d4SSatish Balay   t[i++] = 0.6250;
3019371c9d4SSatish Balay   y[i]   = 60.8000;
3029371c9d4SSatish Balay   t[i++] = 0.7500;
3039371c9d4SSatish Balay   y[i]   = 55.5000;
3049371c9d4SSatish Balay   t[i++] = 0.8750;
3059371c9d4SSatish Balay   y[i]   = 50.3000;
3069371c9d4SSatish Balay   t[i++] = 1.0000;
3079371c9d4SSatish Balay   y[i]   = 41.0000;
3089371c9d4SSatish Balay   t[i++] = 1.2500;
3099371c9d4SSatish Balay   y[i]   = 29.4000;
3109371c9d4SSatish Balay   t[i++] = 1.7500;
3119371c9d4SSatish Balay   y[i]   = 20.4000;
3129371c9d4SSatish Balay   t[i++] = 2.2500;
3139371c9d4SSatish Balay   y[i]   = 29.3625;
3149371c9d4SSatish Balay   t[i++] = 1.7500;
3159371c9d4SSatish Balay   y[i]   = 21.1500;
3169371c9d4SSatish Balay   t[i++] = 2.2500;
3179371c9d4SSatish Balay   y[i]   = 16.7625;
3189371c9d4SSatish Balay   t[i++] = 2.7500;
3199371c9d4SSatish Balay   y[i]   = 13.2000;
3209371c9d4SSatish Balay   t[i++] = 3.2500;
3219371c9d4SSatish Balay   y[i]   = 10.8750;
3229371c9d4SSatish Balay   t[i++] = 3.7500;
3239371c9d4SSatish Balay   y[i]   = 8.1750;
3249371c9d4SSatish Balay   t[i++] = 4.2500;
3259371c9d4SSatish Balay   y[i]   = 7.3500;
3269371c9d4SSatish Balay   t[i++] = 4.7500;
3279371c9d4SSatish Balay   y[i]   = 5.9625;
3289371c9d4SSatish Balay   t[i++] = 5.2500;
3299371c9d4SSatish Balay   y[i]   = 5.6250;
3309371c9d4SSatish Balay   t[i++] = 5.7500;
3319371c9d4SSatish Balay   y[i]   = 81.5000;
3329371c9d4SSatish Balay   t[i++] = .5000;
3339371c9d4SSatish Balay   y[i]   = 62.4000;
3349371c9d4SSatish Balay   t[i++] = .7500;
3359371c9d4SSatish Balay   y[i]   = 32.5000;
3369371c9d4SSatish Balay   t[i++] = 1.5000;
3379371c9d4SSatish Balay   y[i]   = 12.4100;
3389371c9d4SSatish Balay   t[i++] = 3.0000;
3399371c9d4SSatish Balay   y[i]   = 13.1200;
3409371c9d4SSatish Balay   t[i++] = 3.0000;
3419371c9d4SSatish Balay   y[i]   = 15.5600;
3429371c9d4SSatish Balay   t[i++] = 3.0000;
3439371c9d4SSatish Balay   y[i]   = 5.6300;
3449371c9d4SSatish Balay   t[i++] = 6.0000;
3459371c9d4SSatish Balay   y[i]   = 78.0000;
3469371c9d4SSatish Balay   t[i++] = .5000;
3479371c9d4SSatish Balay   y[i]   = 59.9000;
3489371c9d4SSatish Balay   t[i++] = .7500;
3499371c9d4SSatish Balay   y[i]   = 33.2000;
3509371c9d4SSatish Balay   t[i++] = 1.5000;
3519371c9d4SSatish Balay   y[i]   = 13.8400;
3529371c9d4SSatish Balay   t[i++] = 3.0000;
3539371c9d4SSatish Balay   y[i]   = 12.7500;
3549371c9d4SSatish Balay   t[i++] = 3.0000;
3559371c9d4SSatish Balay   y[i]   = 14.6200;
3569371c9d4SSatish Balay   t[i++] = 3.0000;
3579371c9d4SSatish Balay   y[i]   = 3.9400;
3589371c9d4SSatish Balay   t[i++] = 6.0000;
3599371c9d4SSatish Balay   y[i]   = 76.8000;
3609371c9d4SSatish Balay   t[i++] = .5000;
3619371c9d4SSatish Balay   y[i]   = 61.0000;
3629371c9d4SSatish Balay   t[i++] = .7500;
3639371c9d4SSatish Balay   y[i]   = 32.9000;
3649371c9d4SSatish Balay   t[i++] = 1.5000;
3659371c9d4SSatish Balay   y[i]   = 13.8700;
3669371c9d4SSatish Balay   t[i++] = 3.0000;
3679371c9d4SSatish Balay   y[i]   = 11.8100;
3689371c9d4SSatish Balay   t[i++] = 3.0000;
3699371c9d4SSatish Balay   y[i]   = 13.3100;
3709371c9d4SSatish Balay   t[i++] = 3.0000;
3719371c9d4SSatish Balay   y[i]   = 5.4400;
3729371c9d4SSatish Balay   t[i++] = 6.0000;
3739371c9d4SSatish Balay   y[i]   = 78.0000;
3749371c9d4SSatish Balay   t[i++] = .5000;
3759371c9d4SSatish Balay   y[i]   = 63.5000;
3769371c9d4SSatish Balay   t[i++] = .7500;
3779371c9d4SSatish Balay   y[i]   = 33.8000;
3789371c9d4SSatish Balay   t[i++] = 1.5000;
3799371c9d4SSatish Balay   y[i]   = 12.5600;
3809371c9d4SSatish Balay   t[i++] = 3.0000;
3819371c9d4SSatish Balay   y[i]   = 5.6300;
3829371c9d4SSatish Balay   t[i++] = 6.0000;
3839371c9d4SSatish Balay   y[i]   = 12.7500;
3849371c9d4SSatish Balay   t[i++] = 3.0000;
3859371c9d4SSatish Balay   y[i]   = 13.1200;
3869371c9d4SSatish Balay   t[i++] = 3.0000;
3879371c9d4SSatish Balay   y[i]   = 5.4400;
3889371c9d4SSatish Balay   t[i++] = 6.0000;
3899371c9d4SSatish Balay   y[i]   = 76.8000;
3909371c9d4SSatish Balay   t[i++] = .5000;
3919371c9d4SSatish Balay   y[i]   = 60.0000;
3929371c9d4SSatish Balay   t[i++] = .7500;
3939371c9d4SSatish Balay   y[i]   = 47.8000;
3949371c9d4SSatish Balay   t[i++] = 1.0000;
3959371c9d4SSatish Balay   y[i]   = 32.0000;
3969371c9d4SSatish Balay   t[i++] = 1.5000;
3979371c9d4SSatish Balay   y[i]   = 22.2000;
3989371c9d4SSatish Balay   t[i++] = 2.0000;
3999371c9d4SSatish Balay   y[i]   = 22.5700;
4009371c9d4SSatish Balay   t[i++] = 2.0000;
4019371c9d4SSatish Balay   y[i]   = 18.8200;
4029371c9d4SSatish Balay   t[i++] = 2.5000;
4039371c9d4SSatish Balay   y[i]   = 13.9500;
4049371c9d4SSatish Balay   t[i++] = 3.0000;
4059371c9d4SSatish Balay   y[i]   = 11.2500;
4069371c9d4SSatish Balay   t[i++] = 4.0000;
4079371c9d4SSatish Balay   y[i]   = 9.0000;
4089371c9d4SSatish Balay   t[i++] = 5.0000;
4099371c9d4SSatish Balay   y[i]   = 6.6700;
4109371c9d4SSatish Balay   t[i++] = 6.0000;
4119371c9d4SSatish Balay   y[i]   = 75.8000;
4129371c9d4SSatish Balay   t[i++] = .5000;
4139371c9d4SSatish Balay   y[i]   = 62.0000;
4149371c9d4SSatish Balay   t[i++] = .7500;
4159371c9d4SSatish Balay   y[i]   = 48.8000;
4169371c9d4SSatish Balay   t[i++] = 1.0000;
4179371c9d4SSatish Balay   y[i]   = 35.2000;
4189371c9d4SSatish Balay   t[i++] = 1.5000;
4199371c9d4SSatish Balay   y[i]   = 20.0000;
4209371c9d4SSatish Balay   t[i++] = 2.0000;
4219371c9d4SSatish Balay   y[i]   = 20.3200;
4229371c9d4SSatish Balay   t[i++] = 2.0000;
4239371c9d4SSatish Balay   y[i]   = 19.3100;
4249371c9d4SSatish Balay   t[i++] = 2.5000;
4259371c9d4SSatish Balay   y[i]   = 12.7500;
4269371c9d4SSatish Balay   t[i++] = 3.0000;
4279371c9d4SSatish Balay   y[i]   = 10.4200;
4289371c9d4SSatish Balay   t[i++] = 4.0000;
4299371c9d4SSatish Balay   y[i]   = 7.3100;
4309371c9d4SSatish Balay   t[i++] = 5.0000;
4319371c9d4SSatish Balay   y[i]   = 7.4200;
4329371c9d4SSatish Balay   t[i++] = 6.0000;
4339371c9d4SSatish Balay   y[i]   = 70.5000;
4349371c9d4SSatish Balay   t[i++] = .5000;
4359371c9d4SSatish Balay   y[i]   = 59.5000;
4369371c9d4SSatish Balay   t[i++] = .7500;
4379371c9d4SSatish Balay   y[i]   = 48.5000;
4389371c9d4SSatish Balay   t[i++] = 1.0000;
4399371c9d4SSatish Balay   y[i]   = 35.8000;
4409371c9d4SSatish Balay   t[i++] = 1.5000;
4419371c9d4SSatish Balay   y[i]   = 21.0000;
4429371c9d4SSatish Balay   t[i++] = 2.0000;
4439371c9d4SSatish Balay   y[i]   = 21.6700;
4449371c9d4SSatish Balay   t[i++] = 2.0000;
4459371c9d4SSatish Balay   y[i]   = 21.0000;
4469371c9d4SSatish Balay   t[i++] = 2.5000;
4479371c9d4SSatish Balay   y[i]   = 15.6400;
4489371c9d4SSatish Balay   t[i++] = 3.0000;
4499371c9d4SSatish Balay   y[i]   = 8.1700;
4509371c9d4SSatish Balay   t[i++] = 4.0000;
4519371c9d4SSatish Balay   y[i]   = 8.5500;
4529371c9d4SSatish Balay   t[i++] = 5.0000;
4539371c9d4SSatish Balay   y[i]   = 10.1200;
4549371c9d4SSatish Balay   t[i++] = 6.0000;
4559371c9d4SSatish Balay   y[i]   = 78.0000;
4569371c9d4SSatish Balay   t[i++] = .5000;
4579371c9d4SSatish Balay   y[i]   = 66.0000;
4589371c9d4SSatish Balay   t[i++] = .6250;
4599371c9d4SSatish Balay   y[i]   = 62.0000;
4609371c9d4SSatish Balay   t[i++] = .7500;
4619371c9d4SSatish Balay   y[i]   = 58.0000;
4629371c9d4SSatish Balay   t[i++] = .8750;
4639371c9d4SSatish Balay   y[i]   = 47.7000;
4649371c9d4SSatish Balay   t[i++] = 1.0000;
4659371c9d4SSatish Balay   y[i]   = 37.8000;
4669371c9d4SSatish Balay   t[i++] = 1.2500;
4679371c9d4SSatish Balay   y[i]   = 20.2000;
4689371c9d4SSatish Balay   t[i++] = 2.2500;
4699371c9d4SSatish Balay   y[i]   = 21.0700;
4709371c9d4SSatish Balay   t[i++] = 2.2500;
4719371c9d4SSatish Balay   y[i]   = 13.8700;
4729371c9d4SSatish Balay   t[i++] = 2.7500;
4739371c9d4SSatish Balay   y[i]   = 9.6700;
4749371c9d4SSatish Balay   t[i++] = 3.2500;
4759371c9d4SSatish Balay   y[i]   = 7.7600;
4769371c9d4SSatish Balay   t[i++] = 3.7500;
4779371c9d4SSatish Balay   y[i]   = 5.4400;
4789371c9d4SSatish Balay   t[i++] = 4.2500;
4799371c9d4SSatish Balay   y[i]   = 4.8700;
4809371c9d4SSatish Balay   t[i++] = 4.7500;
4819371c9d4SSatish Balay   y[i]   = 4.0100;
4829371c9d4SSatish Balay   t[i++] = 5.2500;
4839371c9d4SSatish Balay   y[i]   = 3.7500;
4849371c9d4SSatish Balay   t[i++] = 5.7500;
4859371c9d4SSatish Balay   y[i]   = 24.1900;
4869371c9d4SSatish Balay   t[i++] = 3.0000;
4879371c9d4SSatish Balay   y[i]   = 25.7600;
4889371c9d4SSatish Balay   t[i++] = 3.0000;
4899371c9d4SSatish Balay   y[i]   = 18.0700;
4909371c9d4SSatish Balay   t[i++] = 3.0000;
4919371c9d4SSatish Balay   y[i]   = 11.8100;
4929371c9d4SSatish Balay   t[i++] = 3.0000;
4939371c9d4SSatish Balay   y[i]   = 12.0700;
4949371c9d4SSatish Balay   t[i++] = 3.0000;
4959371c9d4SSatish Balay   y[i]   = 16.1200;
4969371c9d4SSatish Balay   t[i++] = 3.0000;
4979371c9d4SSatish Balay   y[i]   = 70.8000;
4989371c9d4SSatish Balay   t[i++] = .5000;
4999371c9d4SSatish Balay   y[i]   = 54.7000;
5009371c9d4SSatish Balay   t[i++] = .7500;
5019371c9d4SSatish Balay   y[i]   = 48.0000;
5029371c9d4SSatish Balay   t[i++] = 1.0000;
5039371c9d4SSatish Balay   y[i]   = 39.8000;
5049371c9d4SSatish Balay   t[i++] = 1.5000;
5059371c9d4SSatish Balay   y[i]   = 29.8000;
5069371c9d4SSatish Balay   t[i++] = 2.0000;
5079371c9d4SSatish Balay   y[i]   = 23.7000;
5089371c9d4SSatish Balay   t[i++] = 2.5000;
5099371c9d4SSatish Balay   y[i]   = 29.6200;
5109371c9d4SSatish Balay   t[i++] = 2.0000;
5119371c9d4SSatish Balay   y[i]   = 23.8100;
5129371c9d4SSatish Balay   t[i++] = 2.5000;
5139371c9d4SSatish Balay   y[i]   = 17.7000;
5149371c9d4SSatish Balay   t[i++] = 3.0000;
5159371c9d4SSatish Balay   y[i]   = 11.5500;
5169371c9d4SSatish Balay   t[i++] = 4.0000;
5179371c9d4SSatish Balay   y[i]   = 12.0700;
5189371c9d4SSatish Balay   t[i++] = 5.0000;
5199371c9d4SSatish Balay   y[i]   = 8.7400;
5209371c9d4SSatish Balay   t[i++] = 6.0000;
5219371c9d4SSatish Balay   y[i]   = 80.7000;
5229371c9d4SSatish Balay   t[i++] = .5000;
5239371c9d4SSatish Balay   y[i]   = 61.3000;
5249371c9d4SSatish Balay   t[i++] = .7500;
5259371c9d4SSatish Balay   y[i]   = 47.5000;
5269371c9d4SSatish Balay   t[i++] = 1.0000;
5279371c9d4SSatish Balay   y[i]   = 29.0000;
5289371c9d4SSatish Balay   t[i++] = 1.5000;
5299371c9d4SSatish Balay   y[i]   = 24.0000;
5309371c9d4SSatish Balay   t[i++] = 2.0000;
5319371c9d4SSatish Balay   y[i]   = 17.7000;
5329371c9d4SSatish Balay   t[i++] = 2.5000;
5339371c9d4SSatish Balay   y[i]   = 24.5600;
5349371c9d4SSatish Balay   t[i++] = 2.0000;
5359371c9d4SSatish Balay   y[i]   = 18.6700;
5369371c9d4SSatish Balay   t[i++] = 2.5000;
5379371c9d4SSatish Balay   y[i]   = 16.2400;
5389371c9d4SSatish Balay   t[i++] = 3.0000;
5399371c9d4SSatish Balay   y[i]   = 8.7400;
5409371c9d4SSatish Balay   t[i++] = 4.0000;
5419371c9d4SSatish Balay   y[i]   = 7.8700;
5429371c9d4SSatish Balay   t[i++] = 5.0000;
5439371c9d4SSatish Balay   y[i]   = 8.5100;
5449371c9d4SSatish Balay   t[i++] = 6.0000;
5459371c9d4SSatish Balay   y[i]   = 66.7000;
5469371c9d4SSatish Balay   t[i++] = .5000;
5479371c9d4SSatish Balay   y[i]   = 59.2000;
5489371c9d4SSatish Balay   t[i++] = .7500;
5499371c9d4SSatish Balay   y[i]   = 40.8000;
5509371c9d4SSatish Balay   t[i++] = 1.0000;
5519371c9d4SSatish Balay   y[i]   = 30.7000;
5529371c9d4SSatish Balay   t[i++] = 1.5000;
5539371c9d4SSatish Balay   y[i]   = 25.7000;
5549371c9d4SSatish Balay   t[i++] = 2.0000;
5559371c9d4SSatish Balay   y[i]   = 16.3000;
5569371c9d4SSatish Balay   t[i++] = 2.5000;
5579371c9d4SSatish Balay   y[i]   = 25.9900;
5589371c9d4SSatish Balay   t[i++] = 2.0000;
5599371c9d4SSatish Balay   y[i]   = 16.9500;
5609371c9d4SSatish Balay   t[i++] = 2.5000;
5619371c9d4SSatish Balay   y[i]   = 13.3500;
5629371c9d4SSatish Balay   t[i++] = 3.0000;
5639371c9d4SSatish Balay   y[i]   = 8.6200;
5649371c9d4SSatish Balay   t[i++] = 4.0000;
5659371c9d4SSatish Balay   y[i]   = 7.2000;
5669371c9d4SSatish Balay   t[i++] = 5.0000;
5679371c9d4SSatish Balay   y[i]   = 6.6400;
5689371c9d4SSatish Balay   t[i++] = 6.0000;
5699371c9d4SSatish Balay   y[i]   = 13.6900;
5709371c9d4SSatish Balay   t[i++] = 3.0000;
5719371c9d4SSatish Balay   y[i]   = 81.0000;
5729371c9d4SSatish Balay   t[i++] = .5000;
5739371c9d4SSatish Balay   y[i]   = 64.5000;
5749371c9d4SSatish Balay   t[i++] = .7500;
5759371c9d4SSatish Balay   y[i]   = 35.5000;
5769371c9d4SSatish Balay   t[i++] = 1.5000;
5779371c9d4SSatish Balay   y[i]   = 13.3100;
5789371c9d4SSatish Balay   t[i++] = 3.0000;
5799371c9d4SSatish Balay   y[i]   = 4.8700;
5809371c9d4SSatish Balay   t[i++] = 6.0000;
5819371c9d4SSatish Balay   y[i]   = 12.9400;
5829371c9d4SSatish Balay   t[i++] = 3.0000;
5839371c9d4SSatish Balay   y[i]   = 5.0600;
5849371c9d4SSatish Balay   t[i++] = 6.0000;
5859371c9d4SSatish Balay   y[i]   = 15.1900;
5869371c9d4SSatish Balay   t[i++] = 3.0000;
5879371c9d4SSatish Balay   y[i]   = 14.6200;
5889371c9d4SSatish Balay   t[i++] = 3.0000;
5899371c9d4SSatish Balay   y[i]   = 15.6400;
5909371c9d4SSatish Balay   t[i++] = 3.0000;
5919371c9d4SSatish Balay   y[i]   = 25.5000;
5929371c9d4SSatish Balay   t[i++] = 1.7500;
5939371c9d4SSatish Balay   y[i]   = 25.9500;
5949371c9d4SSatish Balay   t[i++] = 1.7500;
5959371c9d4SSatish Balay   y[i]   = 81.7000;
5969371c9d4SSatish Balay   t[i++] = .5000;
5979371c9d4SSatish Balay   y[i]   = 61.6000;
5989371c9d4SSatish Balay   t[i++] = .7500;
5999371c9d4SSatish Balay   y[i]   = 29.8000;
6009371c9d4SSatish Balay   t[i++] = 1.7500;
6019371c9d4SSatish Balay   y[i]   = 29.8100;
6029371c9d4SSatish Balay   t[i++] = 1.7500;
6039371c9d4SSatish Balay   y[i]   = 17.1700;
6049371c9d4SSatish Balay   t[i++] = 2.7500;
6059371c9d4SSatish Balay   y[i]   = 10.3900;
6069371c9d4SSatish Balay   t[i++] = 3.7500;
6079371c9d4SSatish Balay   y[i]   = 28.4000;
6089371c9d4SSatish Balay   t[i++] = 1.7500;
6099371c9d4SSatish Balay   y[i]   = 28.6900;
6109371c9d4SSatish Balay   t[i++] = 1.7500;
6119371c9d4SSatish Balay   y[i]   = 81.3000;
6129371c9d4SSatish Balay   t[i++] = .5000;
6139371c9d4SSatish Balay   y[i]   = 60.9000;
6149371c9d4SSatish Balay   t[i++] = .7500;
6159371c9d4SSatish Balay   y[i]   = 16.6500;
6169371c9d4SSatish Balay   t[i++] = 2.7500;
6179371c9d4SSatish Balay   y[i]   = 10.0500;
6189371c9d4SSatish Balay   t[i++] = 3.7500;
6199371c9d4SSatish Balay   y[i]   = 28.9000;
6209371c9d4SSatish Balay   t[i++] = 1.7500;
6219371c9d4SSatish Balay   y[i]   = 28.9500;
6229371c9d4SSatish Balay   t[i++] = 1.7500;
6233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
624c4762a1bSJed Brown }
625c4762a1bSJed Brown 
626c4762a1bSJed Brown /*TEST
627c4762a1bSJed Brown 
628c4762a1bSJed Brown      build:
629c4762a1bSJed Brown        requires: !complex
630c4762a1bSJed Brown 
631c4762a1bSJed Brown      test:
63210978b7dSBarry Smith        args: -tao_monitor_short -tao_max_it 100 -tao_type pounders -tao_pounders_delta 0.05 -tao_gatol 1.e-5
633c4762a1bSJed Brown        requires: !single
634c4762a1bSJed Brown        TODO: produces different output for many different systems
635c4762a1bSJed Brown 
636c4762a1bSJed Brown      test:
637c4762a1bSJed Brown        suffix: 2
63810978b7dSBarry Smith        args: -tao_monitor_short -tao_max_it 100 -wtype 1 -tao_type pounders -tao_pounders_delta 0.05 -tao_gatol 1.e-5
639c4762a1bSJed Brown        requires: !single
640c4762a1bSJed Brown        TODO: produces different output for many different systems
641c4762a1bSJed Brown 
642c4762a1bSJed Brown      test:
643c4762a1bSJed Brown        suffix: 3
64410978b7dSBarry Smith        args: -tao_monitor_short -tao_max_it 100 -wtype 2 -tao_type pounders -tao_pounders_delta 0.05 -tao_gatol 1.e-5
645c4762a1bSJed Brown        requires: !single
646c4762a1bSJed Brown        TODO: produces different output for many different systems
647c4762a1bSJed Brown 
648c4762a1bSJed Brown      test:
649c4762a1bSJed Brown        suffix: 4
65010978b7dSBarry Smith        args: -tao_monitor_short -tao_max_it 100 -tao_type pounders -tao_pounders_delta 0.05 -pounders_subsolver_tao_type blmvm -tao_gatol 1.e-5
651c4762a1bSJed Brown        requires: !single
652c4762a1bSJed Brown        TODO: produces different output for many different systems
653c4762a1bSJed Brown 
654c4762a1bSJed Brown  TEST*/
655