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 /*--------------------------------------------------------------------*/ 46d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 47d71ae5a4SJacob Faibussowitsch { 48c4762a1bSJed Brown Vec x, f; /* solution, function */ 49c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 50c4762a1bSJed Brown Tao tao; /* Tao solver context */ 51c4762a1bSJed Brown PetscInt i; /* iteration information */ 52c4762a1bSJed Brown PetscReal hist[100], resid[100]; 53c4762a1bSJed Brown PetscInt lits[100]; 54c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 55c4762a1bSJed Brown 56327415f7SBarry Smith PetscFunctionBeginUser; 57c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 58c4762a1bSJed Brown /* Allocate vectors */ 599566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(MPI_COMM_SELF, NPARAMETERS, &x)); 609566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(MPI_COMM_SELF, NOBSERVATIONS, &f)); 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* Create the Jacobian matrix. */ 639566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(MPI_COMM_SELF, NOBSERVATIONS, NPARAMETERS, NULL, &J)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown for (i = 0; i < NOBSERVATIONS; i++) user.idm[i] = i; 66c4762a1bSJed Brown 67c4762a1bSJed Brown for (i = 0; i < NPARAMETERS; i++) user.idn[i] = i; 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 709566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF, &tao)); 719566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOPOUNDERS)); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* Set the function and Jacobian routines. */ 749566063dSJacob Faibussowitsch PetscCall(InitializeData(&user)); 759566063dSJacob Faibussowitsch PetscCall(FormStartingPoint(x)); 769566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 779566063dSJacob Faibussowitsch PetscCall(TaoSetResidualRoutine(tao, f, EvaluateFunction, (void *)&user)); 789566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianResidualRoutine(tao, J, J, EvaluateJacobian, (void *)&user)); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* Check for any TAO command line arguments */ 819566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(TaoSetConvergenceHistory(tao, hist, resid, 0, lits, 100, PETSC_TRUE)); 84c4762a1bSJed Brown /* Perform the Solve */ 859566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* View the vector; then destroy it. */ 889566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_SELF)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Free TAO data structures */ 919566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* Free PETSc data structures */ 949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&f)); 969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 97c4762a1bSJed Brown 989566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 99b122ec5aSJacob Faibussowitsch return 0; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown /*--------------------------------------------------------------------*/ 103d71ae5a4SJacob Faibussowitsch PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr) 104d71ae5a4SJacob Faibussowitsch { 105c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 106c4762a1bSJed Brown PetscInt i; 107c4762a1bSJed Brown const PetscReal *x; 108c4762a1bSJed Brown PetscReal *y = user->y, *f, *t = user->t; 109c4762a1bSJed Brown 110c4762a1bSJed Brown PetscFunctionBegin; 1119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 1129566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 113c4762a1bSJed Brown 114ad540459SPierre Jolivet for (i = 0; i < NOBSERVATIONS; i++) f[i] = y[i] - PetscExpScalar(-x[0] * t[i]) / (x[1] + x[2] * t[i]); 1159566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 1173ba16761SJacob Faibussowitsch PetscCall(PetscLogFlops(6 * NOBSERVATIONS)); 1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121c4762a1bSJed Brown /*------------------------------------------------------------*/ 122c4762a1bSJed Brown /* J[i][j] = df[i]/dt[j] */ 123d71ae5a4SJacob Faibussowitsch PetscErrorCode EvaluateJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr) 124d71ae5a4SJacob Faibussowitsch { 125c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 126c4762a1bSJed Brown PetscInt i; 127c4762a1bSJed Brown const PetscReal *x; 128c4762a1bSJed Brown PetscReal *t = user->t; 129c4762a1bSJed Brown PetscReal base; 130c4762a1bSJed Brown 131c4762a1bSJed Brown PetscFunctionBegin; 1329566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 133c4762a1bSJed Brown for (i = 0; i < NOBSERVATIONS; i++) { 134c4762a1bSJed Brown base = PetscExpScalar(-x[0] * t[i]) / (x[1] + x[2] * t[i]); 135c4762a1bSJed Brown 136c4762a1bSJed Brown user->j[i][0] = t[i] * base; 137c4762a1bSJed Brown user->j[i][1] = base / (x[1] + x[2] * t[i]); 138c4762a1bSJed Brown user->j[i][2] = base * t[i] / (x[1] + x[2] * t[i]); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* Assemble the matrix */ 1429566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, NOBSERVATIONS, user->idm, NPARAMETERS, user->idn, (PetscReal *)user->j, INSERT_VALUES)); 1439566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1449566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 145c4762a1bSJed Brown 1469566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1473ba16761SJacob Faibussowitsch PetscCall(PetscLogFlops(NOBSERVATIONS * 13)); 1483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* ------------------------------------------------------------ */ 152d71ae5a4SJacob Faibussowitsch PetscErrorCode FormStartingPoint(Vec X) 153d71ae5a4SJacob Faibussowitsch { 154c4762a1bSJed Brown PetscReal *x; 155c4762a1bSJed Brown 156c4762a1bSJed Brown PetscFunctionBegin; 1579566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 158c4762a1bSJed Brown x[0] = 0.15; 159c4762a1bSJed Brown x[1] = 0.008; 160c4762a1bSJed Brown x[2] = 0.010; 1619566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* ---------------------------------------------------------------------- */ 166d71ae5a4SJacob Faibussowitsch PetscErrorCode InitializeData(AppCtx *user) 167d71ae5a4SJacob Faibussowitsch { 168c4762a1bSJed Brown PetscReal *t = user->t, *y = user->y; 169c4762a1bSJed Brown PetscInt i = 0; 170c4762a1bSJed Brown 171c4762a1bSJed Brown PetscFunctionBegin; 1729371c9d4SSatish Balay y[i] = 92.9000; 1739371c9d4SSatish Balay t[i++] = 0.5000; 1749371c9d4SSatish Balay y[i] = 78.7000; 1759371c9d4SSatish Balay t[i++] = 0.6250; 1769371c9d4SSatish Balay y[i] = 64.2000; 1779371c9d4SSatish Balay t[i++] = 0.7500; 1789371c9d4SSatish Balay y[i] = 64.9000; 1799371c9d4SSatish Balay t[i++] = 0.8750; 1809371c9d4SSatish Balay y[i] = 57.1000; 1819371c9d4SSatish Balay t[i++] = 1.0000; 1829371c9d4SSatish Balay y[i] = 43.3000; 1839371c9d4SSatish Balay t[i++] = 1.2500; 1849371c9d4SSatish Balay y[i] = 31.1000; 1859371c9d4SSatish Balay t[i++] = 1.7500; 1869371c9d4SSatish Balay y[i] = 23.6000; 1879371c9d4SSatish Balay t[i++] = 2.2500; 1889371c9d4SSatish Balay y[i] = 31.0500; 1899371c9d4SSatish Balay t[i++] = 1.7500; 1909371c9d4SSatish Balay y[i] = 23.7750; 1919371c9d4SSatish Balay t[i++] = 2.2500; 1929371c9d4SSatish Balay y[i] = 17.7375; 1939371c9d4SSatish Balay t[i++] = 2.7500; 1949371c9d4SSatish Balay y[i] = 13.8000; 1959371c9d4SSatish Balay t[i++] = 3.2500; 1969371c9d4SSatish Balay y[i] = 11.5875; 1979371c9d4SSatish Balay t[i++] = 3.7500; 1989371c9d4SSatish Balay y[i] = 9.4125; 1999371c9d4SSatish Balay t[i++] = 4.2500; 2009371c9d4SSatish Balay y[i] = 7.7250; 2019371c9d4SSatish Balay t[i++] = 4.7500; 2029371c9d4SSatish Balay y[i] = 7.3500; 2039371c9d4SSatish Balay t[i++] = 5.2500; 2049371c9d4SSatish Balay y[i] = 8.0250; 2059371c9d4SSatish Balay t[i++] = 5.7500; 2069371c9d4SSatish Balay y[i] = 90.6000; 2079371c9d4SSatish Balay t[i++] = 0.5000; 2089371c9d4SSatish Balay y[i] = 76.9000; 2099371c9d4SSatish Balay t[i++] = 0.6250; 2109371c9d4SSatish Balay y[i] = 71.6000; 2119371c9d4SSatish Balay t[i++] = 0.7500; 2129371c9d4SSatish Balay y[i] = 63.6000; 2139371c9d4SSatish Balay t[i++] = 0.8750; 2149371c9d4SSatish Balay y[i] = 54.0000; 2159371c9d4SSatish Balay t[i++] = 1.0000; 2169371c9d4SSatish Balay y[i] = 39.2000; 2179371c9d4SSatish Balay t[i++] = 1.2500; 2189371c9d4SSatish Balay y[i] = 29.3000; 2199371c9d4SSatish Balay t[i++] = 1.7500; 2209371c9d4SSatish Balay y[i] = 21.4000; 2219371c9d4SSatish Balay t[i++] = 2.2500; 2229371c9d4SSatish Balay y[i] = 29.1750; 2239371c9d4SSatish Balay t[i++] = 1.7500; 2249371c9d4SSatish Balay y[i] = 22.1250; 2259371c9d4SSatish Balay t[i++] = 2.2500; 2269371c9d4SSatish Balay y[i] = 17.5125; 2279371c9d4SSatish Balay t[i++] = 2.7500; 2289371c9d4SSatish Balay y[i] = 14.2500; 2299371c9d4SSatish Balay t[i++] = 3.2500; 2309371c9d4SSatish Balay y[i] = 9.4500; 2319371c9d4SSatish Balay t[i++] = 3.7500; 2329371c9d4SSatish Balay y[i] = 9.1500; 2339371c9d4SSatish Balay t[i++] = 4.2500; 2349371c9d4SSatish Balay y[i] = 7.9125; 2359371c9d4SSatish Balay t[i++] = 4.7500; 2369371c9d4SSatish Balay y[i] = 8.4750; 2379371c9d4SSatish Balay t[i++] = 5.2500; 2389371c9d4SSatish Balay y[i] = 6.1125; 2399371c9d4SSatish Balay t[i++] = 5.7500; 2409371c9d4SSatish Balay y[i] = 80.0000; 2419371c9d4SSatish Balay t[i++] = 0.5000; 2429371c9d4SSatish Balay y[i] = 79.0000; 2439371c9d4SSatish Balay t[i++] = 0.6250; 2449371c9d4SSatish Balay y[i] = 63.8000; 2459371c9d4SSatish Balay t[i++] = 0.7500; 2469371c9d4SSatish Balay y[i] = 57.2000; 2479371c9d4SSatish Balay t[i++] = 0.8750; 2489371c9d4SSatish Balay y[i] = 53.2000; 2499371c9d4SSatish Balay t[i++] = 1.0000; 2509371c9d4SSatish Balay y[i] = 42.5000; 2519371c9d4SSatish Balay t[i++] = 1.2500; 2529371c9d4SSatish Balay y[i] = 26.8000; 2539371c9d4SSatish Balay t[i++] = 1.7500; 2549371c9d4SSatish Balay y[i] = 20.4000; 2559371c9d4SSatish Balay t[i++] = 2.2500; 2569371c9d4SSatish Balay y[i] = 26.8500; 2579371c9d4SSatish Balay t[i++] = 1.7500; 2589371c9d4SSatish Balay y[i] = 21.0000; 2599371c9d4SSatish Balay t[i++] = 2.2500; 2609371c9d4SSatish Balay y[i] = 16.4625; 2619371c9d4SSatish Balay t[i++] = 2.7500; 2629371c9d4SSatish Balay y[i] = 12.5250; 2639371c9d4SSatish Balay t[i++] = 3.2500; 2649371c9d4SSatish Balay y[i] = 10.5375; 2659371c9d4SSatish Balay t[i++] = 3.7500; 2669371c9d4SSatish Balay y[i] = 8.5875; 2679371c9d4SSatish Balay t[i++] = 4.2500; 2689371c9d4SSatish Balay y[i] = 7.1250; 2699371c9d4SSatish Balay t[i++] = 4.7500; 2709371c9d4SSatish Balay y[i] = 6.1125; 2719371c9d4SSatish Balay t[i++] = 5.2500; 2729371c9d4SSatish Balay y[i] = 5.9625; 2739371c9d4SSatish Balay t[i++] = 5.7500; 2749371c9d4SSatish Balay y[i] = 74.1000; 2759371c9d4SSatish Balay t[i++] = 0.5000; 2769371c9d4SSatish Balay y[i] = 67.3000; 2779371c9d4SSatish Balay t[i++] = 0.6250; 2789371c9d4SSatish Balay y[i] = 60.8000; 2799371c9d4SSatish Balay t[i++] = 0.7500; 2809371c9d4SSatish Balay y[i] = 55.5000; 2819371c9d4SSatish Balay t[i++] = 0.8750; 2829371c9d4SSatish Balay y[i] = 50.3000; 2839371c9d4SSatish Balay t[i++] = 1.0000; 2849371c9d4SSatish Balay y[i] = 41.0000; 2859371c9d4SSatish Balay t[i++] = 1.2500; 2869371c9d4SSatish Balay y[i] = 29.4000; 2879371c9d4SSatish Balay t[i++] = 1.7500; 2889371c9d4SSatish Balay y[i] = 20.4000; 2899371c9d4SSatish Balay t[i++] = 2.2500; 2909371c9d4SSatish Balay y[i] = 29.3625; 2919371c9d4SSatish Balay t[i++] = 1.7500; 2929371c9d4SSatish Balay y[i] = 21.1500; 2939371c9d4SSatish Balay t[i++] = 2.2500; 2949371c9d4SSatish Balay y[i] = 16.7625; 2959371c9d4SSatish Balay t[i++] = 2.7500; 2969371c9d4SSatish Balay y[i] = 13.2000; 2979371c9d4SSatish Balay t[i++] = 3.2500; 2989371c9d4SSatish Balay y[i] = 10.8750; 2999371c9d4SSatish Balay t[i++] = 3.7500; 3009371c9d4SSatish Balay y[i] = 8.1750; 3019371c9d4SSatish Balay t[i++] = 4.2500; 3029371c9d4SSatish Balay y[i] = 7.3500; 3039371c9d4SSatish Balay t[i++] = 4.7500; 3049371c9d4SSatish Balay y[i] = 5.9625; 3059371c9d4SSatish Balay t[i++] = 5.2500; 3069371c9d4SSatish Balay y[i] = 5.6250; 3079371c9d4SSatish Balay t[i++] = 5.7500; 3089371c9d4SSatish Balay y[i] = 81.5000; 3099371c9d4SSatish Balay t[i++] = .5000; 3109371c9d4SSatish Balay y[i] = 62.4000; 3119371c9d4SSatish Balay t[i++] = .7500; 3129371c9d4SSatish Balay y[i] = 32.5000; 3139371c9d4SSatish Balay t[i++] = 1.5000; 3149371c9d4SSatish Balay y[i] = 12.4100; 3159371c9d4SSatish Balay t[i++] = 3.0000; 3169371c9d4SSatish Balay y[i] = 13.1200; 3179371c9d4SSatish Balay t[i++] = 3.0000; 3189371c9d4SSatish Balay y[i] = 15.5600; 3199371c9d4SSatish Balay t[i++] = 3.0000; 3209371c9d4SSatish Balay y[i] = 5.6300; 3219371c9d4SSatish Balay t[i++] = 6.0000; 3229371c9d4SSatish Balay y[i] = 78.0000; 3239371c9d4SSatish Balay t[i++] = .5000; 3249371c9d4SSatish Balay y[i] = 59.9000; 3259371c9d4SSatish Balay t[i++] = .7500; 3269371c9d4SSatish Balay y[i] = 33.2000; 3279371c9d4SSatish Balay t[i++] = 1.5000; 3289371c9d4SSatish Balay y[i] = 13.8400; 3299371c9d4SSatish Balay t[i++] = 3.0000; 3309371c9d4SSatish Balay y[i] = 12.7500; 3319371c9d4SSatish Balay t[i++] = 3.0000; 3329371c9d4SSatish Balay y[i] = 14.6200; 3339371c9d4SSatish Balay t[i++] = 3.0000; 3349371c9d4SSatish Balay y[i] = 3.9400; 3359371c9d4SSatish Balay t[i++] = 6.0000; 3369371c9d4SSatish Balay y[i] = 76.8000; 3379371c9d4SSatish Balay t[i++] = .5000; 3389371c9d4SSatish Balay y[i] = 61.0000; 3399371c9d4SSatish Balay t[i++] = .7500; 3409371c9d4SSatish Balay y[i] = 32.9000; 3419371c9d4SSatish Balay t[i++] = 1.5000; 3429371c9d4SSatish Balay y[i] = 13.8700; 3439371c9d4SSatish Balay t[i++] = 3.0000; 3449371c9d4SSatish Balay y[i] = 11.8100; 3459371c9d4SSatish Balay t[i++] = 3.0000; 3469371c9d4SSatish Balay y[i] = 13.3100; 3479371c9d4SSatish Balay t[i++] = 3.0000; 3489371c9d4SSatish Balay y[i] = 5.4400; 3499371c9d4SSatish Balay t[i++] = 6.0000; 3509371c9d4SSatish Balay y[i] = 78.0000; 3519371c9d4SSatish Balay t[i++] = .5000; 3529371c9d4SSatish Balay y[i] = 63.5000; 3539371c9d4SSatish Balay t[i++] = .7500; 3549371c9d4SSatish Balay y[i] = 33.8000; 3559371c9d4SSatish Balay t[i++] = 1.5000; 3569371c9d4SSatish Balay y[i] = 12.5600; 3579371c9d4SSatish Balay t[i++] = 3.0000; 3589371c9d4SSatish Balay y[i] = 5.6300; 3599371c9d4SSatish Balay t[i++] = 6.0000; 3609371c9d4SSatish Balay y[i] = 12.7500; 3619371c9d4SSatish Balay t[i++] = 3.0000; 3629371c9d4SSatish Balay y[i] = 13.1200; 3639371c9d4SSatish Balay t[i++] = 3.0000; 3649371c9d4SSatish Balay y[i] = 5.4400; 3659371c9d4SSatish Balay t[i++] = 6.0000; 3669371c9d4SSatish Balay y[i] = 76.8000; 3679371c9d4SSatish Balay t[i++] = .5000; 3689371c9d4SSatish Balay y[i] = 60.0000; 3699371c9d4SSatish Balay t[i++] = .7500; 3709371c9d4SSatish Balay y[i] = 47.8000; 3719371c9d4SSatish Balay t[i++] = 1.0000; 3729371c9d4SSatish Balay y[i] = 32.0000; 3739371c9d4SSatish Balay t[i++] = 1.5000; 3749371c9d4SSatish Balay y[i] = 22.2000; 3759371c9d4SSatish Balay t[i++] = 2.0000; 3769371c9d4SSatish Balay y[i] = 22.5700; 3779371c9d4SSatish Balay t[i++] = 2.0000; 3789371c9d4SSatish Balay y[i] = 18.8200; 3799371c9d4SSatish Balay t[i++] = 2.5000; 3809371c9d4SSatish Balay y[i] = 13.9500; 3819371c9d4SSatish Balay t[i++] = 3.0000; 3829371c9d4SSatish Balay y[i] = 11.2500; 3839371c9d4SSatish Balay t[i++] = 4.0000; 3849371c9d4SSatish Balay y[i] = 9.0000; 3859371c9d4SSatish Balay t[i++] = 5.0000; 3869371c9d4SSatish Balay y[i] = 6.6700; 3879371c9d4SSatish Balay t[i++] = 6.0000; 3889371c9d4SSatish Balay y[i] = 75.8000; 3899371c9d4SSatish Balay t[i++] = .5000; 3909371c9d4SSatish Balay y[i] = 62.0000; 3919371c9d4SSatish Balay t[i++] = .7500; 3929371c9d4SSatish Balay y[i] = 48.8000; 3939371c9d4SSatish Balay t[i++] = 1.0000; 3949371c9d4SSatish Balay y[i] = 35.2000; 3959371c9d4SSatish Balay t[i++] = 1.5000; 3969371c9d4SSatish Balay y[i] = 20.0000; 3979371c9d4SSatish Balay t[i++] = 2.0000; 3989371c9d4SSatish Balay y[i] = 20.3200; 3999371c9d4SSatish Balay t[i++] = 2.0000; 4009371c9d4SSatish Balay y[i] = 19.3100; 4019371c9d4SSatish Balay t[i++] = 2.5000; 4029371c9d4SSatish Balay y[i] = 12.7500; 4039371c9d4SSatish Balay t[i++] = 3.0000; 4049371c9d4SSatish Balay y[i] = 10.4200; 4059371c9d4SSatish Balay t[i++] = 4.0000; 4069371c9d4SSatish Balay y[i] = 7.3100; 4079371c9d4SSatish Balay t[i++] = 5.0000; 4089371c9d4SSatish Balay y[i] = 7.4200; 4099371c9d4SSatish Balay t[i++] = 6.0000; 4109371c9d4SSatish Balay y[i] = 70.5000; 4119371c9d4SSatish Balay t[i++] = .5000; 4129371c9d4SSatish Balay y[i] = 59.5000; 4139371c9d4SSatish Balay t[i++] = .7500; 4149371c9d4SSatish Balay y[i] = 48.5000; 4159371c9d4SSatish Balay t[i++] = 1.0000; 4169371c9d4SSatish Balay y[i] = 35.8000; 4179371c9d4SSatish Balay t[i++] = 1.5000; 4189371c9d4SSatish Balay y[i] = 21.0000; 4199371c9d4SSatish Balay t[i++] = 2.0000; 4209371c9d4SSatish Balay y[i] = 21.6700; 4219371c9d4SSatish Balay t[i++] = 2.0000; 4229371c9d4SSatish Balay y[i] = 21.0000; 4239371c9d4SSatish Balay t[i++] = 2.5000; 4249371c9d4SSatish Balay y[i] = 15.6400; 4259371c9d4SSatish Balay t[i++] = 3.0000; 4269371c9d4SSatish Balay y[i] = 8.1700; 4279371c9d4SSatish Balay t[i++] = 4.0000; 4289371c9d4SSatish Balay y[i] = 8.5500; 4299371c9d4SSatish Balay t[i++] = 5.0000; 4309371c9d4SSatish Balay y[i] = 10.1200; 4319371c9d4SSatish Balay t[i++] = 6.0000; 4329371c9d4SSatish Balay y[i] = 78.0000; 4339371c9d4SSatish Balay t[i++] = .5000; 4349371c9d4SSatish Balay y[i] = 66.0000; 4359371c9d4SSatish Balay t[i++] = .6250; 4369371c9d4SSatish Balay y[i] = 62.0000; 4379371c9d4SSatish Balay t[i++] = .7500; 4389371c9d4SSatish Balay y[i] = 58.0000; 4399371c9d4SSatish Balay t[i++] = .8750; 4409371c9d4SSatish Balay y[i] = 47.7000; 4419371c9d4SSatish Balay t[i++] = 1.0000; 4429371c9d4SSatish Balay y[i] = 37.8000; 4439371c9d4SSatish Balay t[i++] = 1.2500; 4449371c9d4SSatish Balay y[i] = 20.2000; 4459371c9d4SSatish Balay t[i++] = 2.2500; 4469371c9d4SSatish Balay y[i] = 21.0700; 4479371c9d4SSatish Balay t[i++] = 2.2500; 4489371c9d4SSatish Balay y[i] = 13.8700; 4499371c9d4SSatish Balay t[i++] = 2.7500; 4509371c9d4SSatish Balay y[i] = 9.6700; 4519371c9d4SSatish Balay t[i++] = 3.2500; 4529371c9d4SSatish Balay y[i] = 7.7600; 4539371c9d4SSatish Balay t[i++] = 3.7500; 4549371c9d4SSatish Balay y[i] = 5.4400; 4559371c9d4SSatish Balay t[i++] = 4.2500; 4569371c9d4SSatish Balay y[i] = 4.8700; 4579371c9d4SSatish Balay t[i++] = 4.7500; 4589371c9d4SSatish Balay y[i] = 4.0100; 4599371c9d4SSatish Balay t[i++] = 5.2500; 4609371c9d4SSatish Balay y[i] = 3.7500; 4619371c9d4SSatish Balay t[i++] = 5.7500; 4629371c9d4SSatish Balay y[i] = 24.1900; 4639371c9d4SSatish Balay t[i++] = 3.0000; 4649371c9d4SSatish Balay y[i] = 25.7600; 4659371c9d4SSatish Balay t[i++] = 3.0000; 4669371c9d4SSatish Balay y[i] = 18.0700; 4679371c9d4SSatish Balay t[i++] = 3.0000; 4689371c9d4SSatish Balay y[i] = 11.8100; 4699371c9d4SSatish Balay t[i++] = 3.0000; 4709371c9d4SSatish Balay y[i] = 12.0700; 4719371c9d4SSatish Balay t[i++] = 3.0000; 4729371c9d4SSatish Balay y[i] = 16.1200; 4739371c9d4SSatish Balay t[i++] = 3.0000; 4749371c9d4SSatish Balay y[i] = 70.8000; 4759371c9d4SSatish Balay t[i++] = .5000; 4769371c9d4SSatish Balay y[i] = 54.7000; 4779371c9d4SSatish Balay t[i++] = .7500; 4789371c9d4SSatish Balay y[i] = 48.0000; 4799371c9d4SSatish Balay t[i++] = 1.0000; 4809371c9d4SSatish Balay y[i] = 39.8000; 4819371c9d4SSatish Balay t[i++] = 1.5000; 4829371c9d4SSatish Balay y[i] = 29.8000; 4839371c9d4SSatish Balay t[i++] = 2.0000; 4849371c9d4SSatish Balay y[i] = 23.7000; 4859371c9d4SSatish Balay t[i++] = 2.5000; 4869371c9d4SSatish Balay y[i] = 29.6200; 4879371c9d4SSatish Balay t[i++] = 2.0000; 4889371c9d4SSatish Balay y[i] = 23.8100; 4899371c9d4SSatish Balay t[i++] = 2.5000; 4909371c9d4SSatish Balay y[i] = 17.7000; 4919371c9d4SSatish Balay t[i++] = 3.0000; 4929371c9d4SSatish Balay y[i] = 11.5500; 4939371c9d4SSatish Balay t[i++] = 4.0000; 4949371c9d4SSatish Balay y[i] = 12.0700; 4959371c9d4SSatish Balay t[i++] = 5.0000; 4969371c9d4SSatish Balay y[i] = 8.7400; 4979371c9d4SSatish Balay t[i++] = 6.0000; 4989371c9d4SSatish Balay y[i] = 80.7000; 4999371c9d4SSatish Balay t[i++] = .5000; 5009371c9d4SSatish Balay y[i] = 61.3000; 5019371c9d4SSatish Balay t[i++] = .7500; 5029371c9d4SSatish Balay y[i] = 47.5000; 5039371c9d4SSatish Balay t[i++] = 1.0000; 5049371c9d4SSatish Balay y[i] = 29.0000; 5059371c9d4SSatish Balay t[i++] = 1.5000; 5069371c9d4SSatish Balay y[i] = 24.0000; 5079371c9d4SSatish Balay t[i++] = 2.0000; 5089371c9d4SSatish Balay y[i] = 17.7000; 5099371c9d4SSatish Balay t[i++] = 2.5000; 5109371c9d4SSatish Balay y[i] = 24.5600; 5119371c9d4SSatish Balay t[i++] = 2.0000; 5129371c9d4SSatish Balay y[i] = 18.6700; 5139371c9d4SSatish Balay t[i++] = 2.5000; 5149371c9d4SSatish Balay y[i] = 16.2400; 5159371c9d4SSatish Balay t[i++] = 3.0000; 5169371c9d4SSatish Balay y[i] = 8.7400; 5179371c9d4SSatish Balay t[i++] = 4.0000; 5189371c9d4SSatish Balay y[i] = 7.8700; 5199371c9d4SSatish Balay t[i++] = 5.0000; 5209371c9d4SSatish Balay y[i] = 8.5100; 5219371c9d4SSatish Balay t[i++] = 6.0000; 5229371c9d4SSatish Balay y[i] = 66.7000; 5239371c9d4SSatish Balay t[i++] = .5000; 5249371c9d4SSatish Balay y[i] = 59.2000; 5259371c9d4SSatish Balay t[i++] = .7500; 5269371c9d4SSatish Balay y[i] = 40.8000; 5279371c9d4SSatish Balay t[i++] = 1.0000; 5289371c9d4SSatish Balay y[i] = 30.7000; 5299371c9d4SSatish Balay t[i++] = 1.5000; 5309371c9d4SSatish Balay y[i] = 25.7000; 5319371c9d4SSatish Balay t[i++] = 2.0000; 5329371c9d4SSatish Balay y[i] = 16.3000; 5339371c9d4SSatish Balay t[i++] = 2.5000; 5349371c9d4SSatish Balay y[i] = 25.9900; 5359371c9d4SSatish Balay t[i++] = 2.0000; 5369371c9d4SSatish Balay y[i] = 16.9500; 5379371c9d4SSatish Balay t[i++] = 2.5000; 5389371c9d4SSatish Balay y[i] = 13.3500; 5399371c9d4SSatish Balay t[i++] = 3.0000; 5409371c9d4SSatish Balay y[i] = 8.6200; 5419371c9d4SSatish Balay t[i++] = 4.0000; 5429371c9d4SSatish Balay y[i] = 7.2000; 5439371c9d4SSatish Balay t[i++] = 5.0000; 5449371c9d4SSatish Balay y[i] = 6.6400; 5459371c9d4SSatish Balay t[i++] = 6.0000; 5469371c9d4SSatish Balay y[i] = 13.6900; 5479371c9d4SSatish Balay t[i++] = 3.0000; 5489371c9d4SSatish Balay y[i] = 81.0000; 5499371c9d4SSatish Balay t[i++] = .5000; 5509371c9d4SSatish Balay y[i] = 64.5000; 5519371c9d4SSatish Balay t[i++] = .7500; 5529371c9d4SSatish Balay y[i] = 35.5000; 5539371c9d4SSatish Balay t[i++] = 1.5000; 5549371c9d4SSatish Balay y[i] = 13.3100; 5559371c9d4SSatish Balay t[i++] = 3.0000; 5569371c9d4SSatish Balay y[i] = 4.8700; 5579371c9d4SSatish Balay t[i++] = 6.0000; 5589371c9d4SSatish Balay y[i] = 12.9400; 5599371c9d4SSatish Balay t[i++] = 3.0000; 5609371c9d4SSatish Balay y[i] = 5.0600; 5619371c9d4SSatish Balay t[i++] = 6.0000; 5629371c9d4SSatish Balay y[i] = 15.1900; 5639371c9d4SSatish Balay t[i++] = 3.0000; 5649371c9d4SSatish Balay y[i] = 14.6200; 5659371c9d4SSatish Balay t[i++] = 3.0000; 5669371c9d4SSatish Balay y[i] = 15.6400; 5679371c9d4SSatish Balay t[i++] = 3.0000; 5689371c9d4SSatish Balay y[i] = 25.5000; 5699371c9d4SSatish Balay t[i++] = 1.7500; 5709371c9d4SSatish Balay y[i] = 25.9500; 5719371c9d4SSatish Balay t[i++] = 1.7500; 5729371c9d4SSatish Balay y[i] = 81.7000; 5739371c9d4SSatish Balay t[i++] = .5000; 5749371c9d4SSatish Balay y[i] = 61.6000; 5759371c9d4SSatish Balay t[i++] = .7500; 5769371c9d4SSatish Balay y[i] = 29.8000; 5779371c9d4SSatish Balay t[i++] = 1.7500; 5789371c9d4SSatish Balay y[i] = 29.8100; 5799371c9d4SSatish Balay t[i++] = 1.7500; 5809371c9d4SSatish Balay y[i] = 17.1700; 5819371c9d4SSatish Balay t[i++] = 2.7500; 5829371c9d4SSatish Balay y[i] = 10.3900; 5839371c9d4SSatish Balay t[i++] = 3.7500; 5849371c9d4SSatish Balay y[i] = 28.4000; 5859371c9d4SSatish Balay t[i++] = 1.7500; 5869371c9d4SSatish Balay y[i] = 28.6900; 5879371c9d4SSatish Balay t[i++] = 1.7500; 5889371c9d4SSatish Balay y[i] = 81.3000; 5899371c9d4SSatish Balay t[i++] = .5000; 5909371c9d4SSatish Balay y[i] = 60.9000; 5919371c9d4SSatish Balay t[i++] = .7500; 5929371c9d4SSatish Balay y[i] = 16.6500; 5939371c9d4SSatish Balay t[i++] = 2.7500; 5949371c9d4SSatish Balay y[i] = 10.0500; 5959371c9d4SSatish Balay t[i++] = 3.7500; 5969371c9d4SSatish Balay y[i] = 28.9000; 5979371c9d4SSatish Balay t[i++] = 1.7500; 5989371c9d4SSatish Balay y[i] = 28.9500; 5999371c9d4SSatish Balay t[i++] = 1.7500; 6003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 601c4762a1bSJed Brown } 602c4762a1bSJed Brown 603c4762a1bSJed Brown /*TEST 604c4762a1bSJed Brown 605c4762a1bSJed Brown build: 606c4762a1bSJed Brown requires: !complex !single 607c4762a1bSJed Brown 608c4762a1bSJed Brown test: 60910978b7dSBarry Smith args: -tao_monitor_short -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-5 610c4762a1bSJed Brown 611c4762a1bSJed Brown test: 612c4762a1bSJed Brown suffix: 2 613*a336c150SZach Atkins args: -tao_monitor_short -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type l2prox -tao_brgn_regularizer_weight 1e-4 -tao_gatol 1.e-5 -tao_brgn_subsolver_tao_monitor_short 614c4762a1bSJed Brown 615c4762a1bSJed Brown test: 616c4762a1bSJed Brown suffix: 3 617*a336c150SZach Atkins args: -tao_monitor_short -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 -tao_brgn_subsolver_tao_monitor_short 618c4762a1bSJed Brown 619cd1c4666STristan Konolige test: 620cd1c4666STristan Konolige suffix: 4 621*a336c150SZach Atkins args: -tao_monitor_short -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type lm -tao_gatol 1.e-5 -tao_brgn_subsolver_tao_type bnls -tao_brgn_subsolver_tao_monitor_short 622cd1c4666STristan Konolige 623c4762a1bSJed Brown TEST*/ 624