xref: /petsc/src/tao/constrained/tutorials/ex1.c (revision c8b0c2b58db46088bcb0f3a2497e96cb9efa794b)
1aad13602SShrirang Abhyankar /* Program usage: mpiexec -n 2 ./ex1 [-help] [all TAO options] */
2aad13602SShrirang Abhyankar 
3aad13602SShrirang Abhyankar /* ----------------------------------------------------------------------
4628da978Sresundermann min f(x) = (x0 - 2)^2 + (x1 - 2)^2 - 2*(x0 + x1)
5aad13602SShrirang Abhyankar s.t.  x0^2 + x1 - 2 = 0
6aad13602SShrirang Abhyankar       0  <= x0^2 - x1 <= 1
7aad13602SShrirang Abhyankar       -1 <= x0, x1 <= 2
8628da978Sresundermann -->
9628da978Sresundermann       g(x)  = 0
10628da978Sresundermann       h(x) >= 0
11628da978Sresundermann       -1 <= x0, x1 <= 2
12628da978Sresundermann where
13628da978Sresundermann       g(x) = x0^2 + x1 - 2
14628da978Sresundermann       h(x) = [x0^2 - x1
15628da978Sresundermann               1 -(x0^2 - x1)]
16aad13602SShrirang Abhyankar ---------------------------------------------------------------------- */
17aad13602SShrirang Abhyankar 
18aad13602SShrirang Abhyankar #include <petsctao.h>
19aad13602SShrirang Abhyankar 
20da81f932SPierre Jolivet static char help[] = "Solves constrained optimization problem using pdipm.\n\
21aad13602SShrirang Abhyankar Input parameters include:\n\
22aad13602SShrirang Abhyankar   -tao_type pdipm    : sets Tao solver\n\
238e58fa1dSresundermann   -no_eq             : removes the equaility constraints from the problem\n\
24*c8b0c2b5SHansol Suh   -no_bound          : removes the bound constraints from the problem\n\
25f560b561SHong Zhang   -init_view         : view initial object setup\n\
26aad13602SShrirang Abhyankar   -snes_fd           : snes with finite difference Jacobian (needed for pdipm)\n\
27628da978Sresundermann   -snes_compare_explicit : compare user Jacobian with finite difference Jacobian \n\
2810978b7dSBarry Smith   -tao_monitor_constraint_norm : convergence monitor with constraint norm \n\
2910978b7dSBarry Smith   -tao_monitor_solution : view exact solution at each iteration\n\
30f560b561SHong Zhang   Note: external package MUMPS is required to run pdipm in parallel. This is designed for a maximum of 2 processors, the code will error if size > 2.\n";
31aad13602SShrirang Abhyankar 
32aad13602SShrirang Abhyankar /*
33628da978Sresundermann    User-defined application context - contains data needed by the application
34aad13602SShrirang Abhyankar */
35aad13602SShrirang Abhyankar typedef struct {
36aad13602SShrirang Abhyankar   PetscInt   n;  /* Global length of x */
37aad13602SShrirang Abhyankar   PetscInt   ne; /* Global number of equality constraints */
38aad13602SShrirang Abhyankar   PetscInt   ni; /* Global number of inequality constraints */
39*c8b0c2b5SHansol Suh   PetscBool  noeqflag, noboundflag, initview;
40aad13602SShrirang Abhyankar   Vec        x, xl, xu;
41aad13602SShrirang Abhyankar   Vec        ce, ci, bl, bu, Xseq;
42aad13602SShrirang Abhyankar   Mat        Ae, Ai, H;
43aad13602SShrirang Abhyankar   VecScatter scat;
44aad13602SShrirang Abhyankar } AppCtx;
45aad13602SShrirang Abhyankar 
46aad13602SShrirang Abhyankar /* -------- User-defined Routines --------- */
47aad13602SShrirang Abhyankar PetscErrorCode InitializeProblem(AppCtx *);
48aad13602SShrirang Abhyankar PetscErrorCode DestroyProblem(AppCtx *);
49aad13602SShrirang Abhyankar PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
5083f04556SHong Zhang PetscErrorCode FormPDIPMHessian(Tao, Vec, Mat, Mat, void *);
51aad13602SShrirang Abhyankar PetscErrorCode FormInequalityConstraints(Tao, Vec, Vec, void *);
52aad13602SShrirang Abhyankar PetscErrorCode FormEqualityConstraints(Tao, Vec, Vec, void *);
53aad13602SShrirang Abhyankar PetscErrorCode FormInequalityJacobian(Tao, Vec, Mat, Mat, void *);
54aad13602SShrirang Abhyankar PetscErrorCode FormEqualityJacobian(Tao, Vec, Mat, Mat, void *);
55aad13602SShrirang Abhyankar 
563ba16761SJacob Faibussowitsch int main(int argc, char **argv)
57d71ae5a4SJacob Faibussowitsch {
58aad13602SShrirang Abhyankar   Tao         tao;
59aad13602SShrirang Abhyankar   KSP         ksp;
60aad13602SShrirang Abhyankar   PC          pc;
61aad13602SShrirang Abhyankar   AppCtx      user; /* application context */
62f560b561SHong Zhang   Vec         x, G, CI, CE;
63f560b561SHong Zhang   PetscMPIInt size;
64661095bbSAlp Dener   TaoType     type;
65f560b561SHong Zhang   PetscReal   f;
669a09327dSAlp Dener   PetscBool   pdipm, mumps;
67aad13602SShrirang Abhyankar 
68327415f7SBarry Smith   PetscFunctionBeginUser;
69c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
709566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
719a09327dSAlp Dener   PetscCheck(size <= 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "More than 2 processors detected. Example written to use max of 2 processors.");
72aad13602SShrirang Abhyankar 
739566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "---- Constrained Problem -----\n"));
749566063dSJacob Faibussowitsch   PetscCall(InitializeProblem(&user)); /* sets up problem, function below */
759a09327dSAlp Dener 
769566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
779a09327dSAlp Dener   PetscCall(TaoSetType(tao, TAOALMM));
789a09327dSAlp Dener   PetscCall(TaoSetSolution(tao, user.x));
79*c8b0c2b5SHansol Suh   if (!user.noboundflag) PetscCall(TaoSetVariableBounds(tao, user.xl, user.xu));
809566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
819a09327dSAlp Dener   PetscCall(TaoSetTolerances(tao, 1.e-4, 0.0, 0.0));
829a09327dSAlp Dener   PetscCall(TaoSetConstraintTolerances(tao, 1.e-4, 0.0));
839a09327dSAlp Dener   PetscCall(TaoSetMaximumFunctionEvaluations(tao, 1e6));
849566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
859a09327dSAlp Dener 
869a09327dSAlp Dener   if (!user.noeqflag) {
879a09327dSAlp Dener     PetscCall(TaoSetEqualityConstraintsRoutine(tao, user.ce, FormEqualityConstraints, (void *)&user));
889a09327dSAlp Dener     PetscCall(TaoSetJacobianEqualityRoutine(tao, user.Ae, user.Ae, FormEqualityJacobian, (void *)&user));
899a09327dSAlp Dener   }
909a09327dSAlp Dener   PetscCall(TaoSetInequalityConstraintsRoutine(tao, user.ci, FormInequalityConstraints, (void *)&user));
919a09327dSAlp Dener   PetscCall(TaoSetJacobianInequalityRoutine(tao, user.Ai, user.Ai, FormInequalityJacobian, (void *)&user));
929a09327dSAlp Dener 
939566063dSJacob Faibussowitsch   PetscCall(TaoGetType(tao, &type));
949566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOPDIPM, &pdipm));
95661095bbSAlp Dener   if (pdipm) {
969a09327dSAlp Dener     /*
979a09327dSAlp Dener       PDIPM produces an indefinite KKT matrix with some zeros along the diagonal
989a09327dSAlp Dener       Inverting this indefinite matrix requires PETSc to be configured with MUMPS
999a09327dSAlp Dener     */
1009a09327dSAlp Dener     PetscCall(PetscHasExternalPackage("mumps", &mumps));
1019a09327dSAlp Dener     PetscCheck(mumps, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "TAOPDIPM requires PETSc to be configured with MUMPS (--download-mumps)");
1029a09327dSAlp Dener     PetscCall(TaoGetKSP(tao, &ksp));
1039a09327dSAlp Dener     PetscCall(KSPGetPC(ksp, &pc));
1049a09327dSAlp Dener     PetscCall(PCFactorSetMatSolverType(pc, MATSOLVERMUMPS));
1059a09327dSAlp Dener     PetscCall(KSPSetType(ksp, KSPPREONLY));
10683f04556SHong Zhang     PetscCall(PCSetType(pc, PCCHOLESKY));
1079a09327dSAlp Dener     PetscCall(KSPSetFromOptions(ksp));
10883f04556SHong Zhang     PetscCall(TaoSetHessian(tao, user.H, user.H, FormPDIPMHessian, (void *)&user));
1099a09327dSAlp Dener   }
1109a09327dSAlp Dener 
1119a09327dSAlp Dener   /* Print out an initial view of the problem */
112f560b561SHong Zhang   if (user.initview) {
1139566063dSJacob Faibussowitsch     PetscCall(TaoSetUp(tao));
1149566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(user.x, &G));
1159566063dSJacob Faibussowitsch     PetscCall(FormFunctionGradient(tao, user.x, &f, G, (void *)&user));
116*c8b0c2b5SHansol Suh     if (pdipm) PetscCall(FormPDIPMHessian(tao, user.x, user.H, user.H, (void *)&user));
1179566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_WORLD));
11863a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial point X:\n"));
1199566063dSJacob Faibussowitsch     PetscCall(VecView(user.x, PETSC_VIEWER_STDOUT_WORLD));
12063a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial objective f(x) = %g\n", (double)f));
12163a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial gradient and Hessian:\n"));
1229566063dSJacob Faibussowitsch     PetscCall(VecView(G, PETSC_VIEWER_STDOUT_WORLD));
123*c8b0c2b5SHansol Suh     if (pdipm) PetscCall(MatView(user.H, PETSC_VIEWER_STDOUT_WORLD));
1249566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&G));
1259566063dSJacob Faibussowitsch     PetscCall(FormInequalityJacobian(tao, user.x, user.Ai, user.Ai, (void *)&user));
1269566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(user.Ai, NULL, &CI));
1279566063dSJacob Faibussowitsch     PetscCall(FormInequalityConstraints(tao, user.x, CI, (void *)&user));
12863a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial inequality constraints and Jacobian:\n"));
1299566063dSJacob Faibussowitsch     PetscCall(VecView(CI, PETSC_VIEWER_STDOUT_WORLD));
1309566063dSJacob Faibussowitsch     PetscCall(MatView(user.Ai, PETSC_VIEWER_STDOUT_WORLD));
1319566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&CI));
132f560b561SHong Zhang     if (!user.noeqflag) {
1339566063dSJacob Faibussowitsch       PetscCall(FormEqualityJacobian(tao, user.x, user.Ae, user.Ae, (void *)&user));
1349566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(user.Ae, NULL, &CE));
1359566063dSJacob Faibussowitsch       PetscCall(FormEqualityConstraints(tao, user.x, CE, (void *)&user));
13663a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial equality constraints and Jacobian:\n"));
1379566063dSJacob Faibussowitsch       PetscCall(VecView(CE, PETSC_VIEWER_STDOUT_WORLD));
1389566063dSJacob Faibussowitsch       PetscCall(MatView(user.Ae, PETSC_VIEWER_STDOUT_WORLD));
1399566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&CE));
140f560b561SHong Zhang     }
1419566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
1429566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_WORLD));
143f560b561SHong Zhang   }
144aad13602SShrirang Abhyankar 
1459566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
1469566063dSJacob Faibussowitsch   PetscCall(TaoGetSolution(tao, &x));
1479a09327dSAlp Dener   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Found solution:\n"));
1489566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
149aad13602SShrirang Abhyankar 
150aad13602SShrirang Abhyankar   /* Free objects */
1519566063dSJacob Faibussowitsch   PetscCall(DestroyProblem(&user));
1529566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
1539566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1543ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
155aad13602SShrirang Abhyankar }
156aad13602SShrirang Abhyankar 
157d71ae5a4SJacob Faibussowitsch PetscErrorCode InitializeProblem(AppCtx *user)
158d71ae5a4SJacob Faibussowitsch {
159aad13602SShrirang Abhyankar   PetscMPIInt size;
160aad13602SShrirang Abhyankar   PetscMPIInt rank;
161aad13602SShrirang Abhyankar   PetscInt    nloc, neloc, niloc;
162aad13602SShrirang Abhyankar 
163aad13602SShrirang Abhyankar   PetscFunctionBegin;
1649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1668e58fa1dSresundermann   user->noeqflag    = PETSC_FALSE;
167*c8b0c2b5SHansol Suh   user->noboundflag = PETSC_FALSE;
1689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_eq", &user->noeqflag, NULL));
169*c8b0c2b5SHansol Suh   PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_bound", &user->noboundflag, NULL));
170f560b561SHong Zhang   user->initview = PETSC_FALSE;
1719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-init_view", &user->initview, NULL));
172628da978Sresundermann 
173f560b561SHong Zhang   /* Tell user the correct solution, not an error checking */
1749a09327dSAlp Dener   if (!user->noeqflag) {
175*c8b0c2b5SHansol Suh     /* With equality constraint, bound or no bound makes no different to the solution */
1769a09327dSAlp Dener     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Expected solution: f(1, 1) = -2\n"));
1779a09327dSAlp Dener   } else {
178*c8b0c2b5SHansol Suh     if (user->noboundflag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Expected solution (-no_eq, -no_bound): f(2.05655, 3.22938) = -9.05728\n"));
179*c8b0c2b5SHansol Suh     else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Expected solution (-no_eq): f(1.73205, 2) = -7.3923\n"));
1808e58fa1dSresundermann   }
181aad13602SShrirang Abhyankar 
1822d4ee042Sprj-   /* create vector x and set initial values */
183aad13602SShrirang Abhyankar   user->n = 2; /* global length */
184f560b561SHong Zhang   nloc    = (size == 1) ? user->n : 1;
1859566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->x));
1869566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->x, nloc, user->n));
1879566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->x));
1889a09327dSAlp Dener   PetscCall(VecSet(user->x, 0.0));
189aad13602SShrirang Abhyankar 
190aad13602SShrirang Abhyankar   /* create and set lower and upper bound vectors */
1919566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->x, &user->xl));
1929566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->x, &user->xu));
1939566063dSJacob Faibussowitsch   PetscCall(VecSet(user->xl, -1.0));
1949566063dSJacob Faibussowitsch   PetscCall(VecSet(user->xu, 2.0));
195aad13602SShrirang Abhyankar 
196aad13602SShrirang Abhyankar   /* create scater to zero */
1979566063dSJacob Faibussowitsch   PetscCall(VecScatterCreateToZero(user->x, &user->scat, &user->Xseq));
198aad13602SShrirang Abhyankar 
199aad13602SShrirang Abhyankar   user->ne = 1;
200628da978Sresundermann   user->ni = 2;
201aad13602SShrirang Abhyankar   neloc    = (rank == 0) ? user->ne : 0;
202f560b561SHong Zhang   niloc    = (size == 1) ? user->ni : 1;
203aad13602SShrirang Abhyankar 
2048e58fa1dSresundermann   if (!user->noeqflag) {
2059566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_WORLD, &user->ce)); /* a 1x1 vec for equality constraints */
2069566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(user->ce, neloc, user->ne));
2079566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(user->ce));
2089566063dSJacob Faibussowitsch     PetscCall(VecSetUp(user->ce));
2098e58fa1dSresundermann   }
210628da978Sresundermann 
2119566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->ci)); /* a 2x1 vec for inequality constraints */
2129566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->ci, niloc, user->ni));
2139566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->ci));
2149566063dSJacob Faibussowitsch   PetscCall(VecSetUp(user->ci));
215aad13602SShrirang Abhyankar 
216baca6076SPierre Jolivet   /* nexn & nixn matrices for equally and inequalty constraints */
2178e58fa1dSresundermann   if (!user->noeqflag) {
2189566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Ae));
2199566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(user->Ae, neloc, nloc, user->ne, user->n));
2209566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(user->Ae));
2219566063dSJacob Faibussowitsch     PetscCall(MatSetUp(user->Ae));
2228e58fa1dSresundermann   }
2238e58fa1dSresundermann 
2249566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Ai));
2259566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Ai, niloc, nloc, user->ni, user->n));
2269566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Ai));
2279566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user->Ai));
228628da978Sresundermann 
2299566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->H));
2309566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->H, nloc, nloc, user->n, user->n));
2319566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->H));
2329566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user->H));
2333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
234aad13602SShrirang Abhyankar }
235aad13602SShrirang Abhyankar 
236d71ae5a4SJacob Faibussowitsch PetscErrorCode DestroyProblem(AppCtx *user)
237d71ae5a4SJacob Faibussowitsch {
238aad13602SShrirang Abhyankar   PetscFunctionBegin;
23948a46eb9SPierre Jolivet   if (!user->noeqflag) PetscCall(MatDestroy(&user->Ae));
2409566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Ai));
2419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->H));
242aad13602SShrirang Abhyankar 
2439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->x));
24448a46eb9SPierre Jolivet   if (!user->noeqflag) PetscCall(VecDestroy(&user->ce));
2459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ci));
2469566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->xl));
2479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->xu));
2489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->Xseq));
2499566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&user->scat));
2503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
251aad13602SShrirang Abhyankar }
252aad13602SShrirang Abhyankar 
253628da978Sresundermann /* Evaluate
254628da978Sresundermann    f(x) = (x0 - 2)^2 + (x1 - 2)^2 - 2*(x0 + x1)
25583f04556SHong Zhang    G = grad f = [2*(x0 - 2) - 2, 2*(x1 - 2) - 2] = 2*X - 6
256aad13602SShrirang Abhyankar */
257d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ctx)
258d71ae5a4SJacob Faibussowitsch {
259aad13602SShrirang Abhyankar   const PetscScalar *x;
260aad13602SShrirang Abhyankar   MPI_Comm           comm;
261aad13602SShrirang Abhyankar   PetscMPIInt        rank;
262aad13602SShrirang Abhyankar   PetscReal          fin;
263aad13602SShrirang Abhyankar   AppCtx            *user = (AppCtx *)ctx;
264aad13602SShrirang Abhyankar   Vec                Xseq = user->Xseq;
265aad13602SShrirang Abhyankar   VecScatter         scat = user->scat;
266aad13602SShrirang Abhyankar 
267aad13602SShrirang Abhyankar   PetscFunctionBegin;
26883f04556SHong Zhang   /* f = (x0 - 2)^2 + (x1 - 2)^2 - 2*(x0 + x1) */
2699566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
2709566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
271aad13602SShrirang Abhyankar 
2729566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
2739566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
274aad13602SShrirang Abhyankar 
275aad13602SShrirang Abhyankar   fin = 0.0;
276dd400576SPatrick Sanan   if (rank == 0) {
2779566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Xseq, &x));
278aad13602SShrirang Abhyankar     fin = (x[0] - 2.0) * (x[0] - 2.0) + (x[1] - 2.0) * (x[1] - 2.0) - 2.0 * (x[0] + x[1]);
2799566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Xseq, &x));
280aad13602SShrirang Abhyankar   }
281462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&fin, f, 1, MPIU_REAL, MPIU_SUM, comm));
28283f04556SHong Zhang 
28383f04556SHong Zhang   /* G = 2*X - 6 */
28483f04556SHong Zhang   PetscCall(VecSet(G, -6.0));
28583f04556SHong Zhang   PetscCall(VecAXPY(G, 2.0, X));
2863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
287aad13602SShrirang Abhyankar }
288aad13602SShrirang Abhyankar 
28983f04556SHong Zhang /* Evaluate PDIPM Hessian, see Eqn(22) in http://doi.org/10.1049/gtd2.12708
29083f04556SHong Zhang    H = Wxx = fxx        + Jacobian(grad g^T*DE)   - Jacobian(grad h^T*DI)]
29183f04556SHong Zhang            = fxx        + Jacobin([2*x0; 1]de[0]) - Jacobian([2*x0, -2*x0; -1, 1][di[0] di[1]]^T)
29283f04556SHong Zhang            = [2 0; 0 2] + [2*de[0]  0;      0  0] - [2*di[0]-2*di[1], 0; 0, 0]
293628da978Sresundermann            = [ 2*(1+de[0]-di[0]+di[1]), 0;
294628da978Sresundermann                           0,            2]
295628da978Sresundermann */
29683f04556SHong Zhang PetscErrorCode FormPDIPMHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx)
297d71ae5a4SJacob Faibussowitsch {
298aad13602SShrirang Abhyankar   AppCtx            *user = (AppCtx *)ctx;
299aad13602SShrirang Abhyankar   Vec                DE, DI;
300aad13602SShrirang Abhyankar   const PetscScalar *de, *di;
301aad13602SShrirang Abhyankar   PetscInt           zero = 0, one = 1;
302aad13602SShrirang Abhyankar   PetscScalar        two = 2.0;
303aad13602SShrirang Abhyankar   PetscScalar        val = 0.0;
304f560b561SHong Zhang   Vec                Deseq, Diseq;
305f560b561SHong Zhang   VecScatter         Descat, Discat;
306aad13602SShrirang Abhyankar   PetscMPIInt        rank;
307aad13602SShrirang Abhyankar   MPI_Comm           comm;
308aad13602SShrirang Abhyankar 
309aad13602SShrirang Abhyankar   PetscFunctionBegin;
3109566063dSJacob Faibussowitsch   PetscCall(TaoGetDualVariables(tao, &DE, &DI));
311aad13602SShrirang Abhyankar 
3129566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
3139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
314aad13602SShrirang Abhyankar 
3158e58fa1dSresundermann   if (!user->noeqflag) {
3169566063dSJacob Faibussowitsch     PetscCall(VecScatterCreateToZero(DE, &Descat, &Deseq));
3179566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(Descat, DE, Deseq, INSERT_VALUES, SCATTER_FORWARD));
3189566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(Descat, DE, Deseq, INSERT_VALUES, SCATTER_FORWARD));
3198e58fa1dSresundermann   }
3209566063dSJacob Faibussowitsch   PetscCall(VecScatterCreateToZero(DI, &Discat, &Diseq));
3219566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(Discat, DI, Diseq, INSERT_VALUES, SCATTER_FORWARD));
3229566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(Discat, DI, Diseq, INSERT_VALUES, SCATTER_FORWARD));
323aad13602SShrirang Abhyankar 
324dd400576SPatrick Sanan   if (rank == 0) {
3259371c9d4SSatish Balay     if (!user->noeqflag) { PetscCall(VecGetArrayRead(Deseq, &de)); /* places equality constraint dual into array */ }
3269566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Diseq, &di)); /* places inequality constraint dual into array */
327628da978Sresundermann 
3288e58fa1dSresundermann     if (!user->noeqflag) {
329628da978Sresundermann       val = 2.0 * (1 + de[0] - di[0] + di[1]);
3309566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(Deseq, &de));
3319566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(Diseq, &di));
3328e58fa1dSresundermann     } else {
333628da978Sresundermann       val = 2.0 * (1 - di[0] + di[1]);
3348e58fa1dSresundermann     }
3359566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Diseq, &di));
3369566063dSJacob Faibussowitsch     PetscCall(MatSetValues(H, 1, &zero, 1, &zero, &val, INSERT_VALUES));
3379566063dSJacob Faibussowitsch     PetscCall(MatSetValues(H, 1, &one, 1, &one, &two, INSERT_VALUES));
338aad13602SShrirang Abhyankar   }
3399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
3409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
3418e58fa1dSresundermann   if (!user->noeqflag) {
3429566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&Descat));
3439566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&Deseq));
3448e58fa1dSresundermann   }
3459566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&Discat));
3469566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Diseq));
3473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
348aad13602SShrirang Abhyankar }
349aad13602SShrirang Abhyankar 
350628da978Sresundermann /* Evaluate
351628da978Sresundermann    h = [ x0^2 - x1;
352628da978Sresundermann          1 -(x0^2 - x1)]
353628da978Sresundermann */
354d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInequalityConstraints(Tao tao, Vec X, Vec CI, void *ctx)
355d71ae5a4SJacob Faibussowitsch {
356aad13602SShrirang Abhyankar   const PetscScalar *x;
357aad13602SShrirang Abhyankar   PetscScalar        ci;
358aad13602SShrirang Abhyankar   MPI_Comm           comm;
359aad13602SShrirang Abhyankar   PetscMPIInt        rank;
360aad13602SShrirang Abhyankar   AppCtx            *user = (AppCtx *)ctx;
361aad13602SShrirang Abhyankar   Vec                Xseq = user->Xseq;
362aad13602SShrirang Abhyankar   VecScatter         scat = user->scat;
363aad13602SShrirang Abhyankar 
364aad13602SShrirang Abhyankar   PetscFunctionBegin;
3659566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
3669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
367aad13602SShrirang Abhyankar 
3689566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
3699566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
370aad13602SShrirang Abhyankar 
371dd400576SPatrick Sanan   if (rank == 0) {
3729566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Xseq, &x));
373aad13602SShrirang Abhyankar     ci = x[0] * x[0] - x[1];
3749566063dSJacob Faibussowitsch     PetscCall(VecSetValue(CI, 0, ci, INSERT_VALUES));
375aad13602SShrirang Abhyankar     ci = -x[0] * x[0] + x[1] + 1.0;
3769566063dSJacob Faibussowitsch     PetscCall(VecSetValue(CI, 1, ci, INSERT_VALUES));
3779566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Xseq, &x));
378aad13602SShrirang Abhyankar   }
3799566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(CI));
3809566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(CI));
3813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
382aad13602SShrirang Abhyankar }
383aad13602SShrirang Abhyankar 
384628da978Sresundermann /* Evaluate
385628da978Sresundermann    g = [ x0^2 + x1 - 2]
386628da978Sresundermann */
387d71ae5a4SJacob Faibussowitsch PetscErrorCode FormEqualityConstraints(Tao tao, Vec X, Vec CE, void *ctx)
388d71ae5a4SJacob Faibussowitsch {
389aad13602SShrirang Abhyankar   const PetscScalar *x;
390aad13602SShrirang Abhyankar   PetscScalar        ce;
391aad13602SShrirang Abhyankar   MPI_Comm           comm;
392aad13602SShrirang Abhyankar   PetscMPIInt        rank;
393aad13602SShrirang Abhyankar   AppCtx            *user = (AppCtx *)ctx;
394aad13602SShrirang Abhyankar   Vec                Xseq = user->Xseq;
395aad13602SShrirang Abhyankar   VecScatter         scat = user->scat;
396aad13602SShrirang Abhyankar 
397aad13602SShrirang Abhyankar   PetscFunctionBegin;
3989566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
3999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
400aad13602SShrirang Abhyankar 
4019566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
4029566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
403aad13602SShrirang Abhyankar 
404dd400576SPatrick Sanan   if (rank == 0) {
4059566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Xseq, &x));
406aad13602SShrirang Abhyankar     ce = x[0] * x[0] + x[1] - 2.0;
4079566063dSJacob Faibussowitsch     PetscCall(VecSetValue(CE, 0, ce, INSERT_VALUES));
4089566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Xseq, &x));
409aad13602SShrirang Abhyankar   }
4109566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(CE));
4119566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(CE));
4123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
413aad13602SShrirang Abhyankar }
414aad13602SShrirang Abhyankar 
415628da978Sresundermann /*
416628da978Sresundermann   grad h = [  2*x0, -1;
417628da978Sresundermann              -2*x0,  1]
418628da978Sresundermann */
419d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInequalityJacobian(Tao tao, Vec X, Mat JI, Mat JIpre, void *ctx)
420d71ae5a4SJacob Faibussowitsch {
421aad13602SShrirang Abhyankar   AppCtx            *user = (AppCtx *)ctx;
422f560b561SHong Zhang   PetscInt           zero = 0, one = 1, cols[2];
423aad13602SShrirang Abhyankar   PetscScalar        vals[2];
424aad13602SShrirang Abhyankar   const PetscScalar *x;
425aad13602SShrirang Abhyankar   Vec                Xseq = user->Xseq;
426aad13602SShrirang Abhyankar   VecScatter         scat = user->scat;
427aad13602SShrirang Abhyankar   MPI_Comm           comm;
428aad13602SShrirang Abhyankar   PetscMPIInt        rank;
429aad13602SShrirang Abhyankar 
430aad13602SShrirang Abhyankar   PetscFunctionBegin;
4319566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
4329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
4339566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
4349566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
435aad13602SShrirang Abhyankar 
4369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xseq, &x));
437c5853193SPierre Jolivet   if (rank == 0) {
4389371c9d4SSatish Balay     cols[0] = 0;
4399371c9d4SSatish Balay     cols[1] = 1;
4409371c9d4SSatish Balay     vals[0] = 2 * x[0];
4419371c9d4SSatish Balay     vals[1] = -1.0;
4429566063dSJacob Faibussowitsch     PetscCall(MatSetValues(JI, 1, &zero, 2, cols, vals, INSERT_VALUES));
4439371c9d4SSatish Balay     vals[0] = -2 * x[0];
4449371c9d4SSatish Balay     vals[1] = 1.0;
4459566063dSJacob Faibussowitsch     PetscCall(MatSetValues(JI, 1, &one, 2, cols, vals, INSERT_VALUES));
446aad13602SShrirang Abhyankar   }
4479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xseq, &x));
4489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(JI, MAT_FINAL_ASSEMBLY));
4499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(JI, MAT_FINAL_ASSEMBLY));
4503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
451aad13602SShrirang Abhyankar }
452aad13602SShrirang Abhyankar 
453628da978Sresundermann /*
45483f04556SHong Zhang   grad g = [2*x0, 1.0]
455628da978Sresundermann */
456d71ae5a4SJacob Faibussowitsch PetscErrorCode FormEqualityJacobian(Tao tao, Vec X, Mat JE, Mat JEpre, void *ctx)
457d71ae5a4SJacob Faibussowitsch {
458f560b561SHong Zhang   PetscInt           zero = 0, cols[2];
459aad13602SShrirang Abhyankar   PetscScalar        vals[2];
460aad13602SShrirang Abhyankar   const PetscScalar *x;
461aad13602SShrirang Abhyankar   PetscMPIInt        rank;
462aad13602SShrirang Abhyankar   MPI_Comm           comm;
463aad13602SShrirang Abhyankar 
464aad13602SShrirang Abhyankar   PetscFunctionBegin;
4659566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
4669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
467aad13602SShrirang Abhyankar 
468dd400576SPatrick Sanan   if (rank == 0) {
4699566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(X, &x));
4709371c9d4SSatish Balay     cols[0] = 0;
4719371c9d4SSatish Balay     cols[1] = 1;
4729371c9d4SSatish Balay     vals[0] = 2 * x[0];
4739371c9d4SSatish Balay     vals[1] = 1.0;
4749566063dSJacob Faibussowitsch     PetscCall(MatSetValues(JE, 1, &zero, 2, cols, vals, INSERT_VALUES));
4759566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(X, &x));
476aad13602SShrirang Abhyankar   }
4779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(JE, MAT_FINAL_ASSEMBLY));
4789566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(JE, MAT_FINAL_ASSEMBLY));
4793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
480aad13602SShrirang Abhyankar }
481aad13602SShrirang Abhyankar 
482aad13602SShrirang Abhyankar /*TEST
483aad13602SShrirang Abhyankar 
484aad13602SShrirang Abhyankar    build:
485f560b561SHong Zhang       requires: !complex !defined(PETSC_USE_CXX)
486aad13602SShrirang Abhyankar 
487aad13602SShrirang Abhyankar    test:
4889a09327dSAlp Dener       args: -tao_converged_reason -tao_gatol 1.e-6 -tao_type pdipm -tao_pdipm_kkt_shift_pd
489910b42cbSStefano Zampini       requires: mumps !single
4909a09327dSAlp Dener       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
491aad13602SShrirang Abhyankar 
492aad13602SShrirang Abhyankar    test:
49383f04556SHong Zhang       suffix: pdipm_2
49483f04556SHong Zhang       requires: mumps
49583f04556SHong Zhang       nsize: 2
49683f04556SHong Zhang       args: -tao_converged_reason -tao_gatol 1.e-6 -tao_type pdipm
49783f04556SHong Zhang       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
49883f04556SHong Zhang 
49983f04556SHong Zhang    test:
500aad13602SShrirang Abhyankar       suffix: 2
5019a09327dSAlp Dener       args: -tao_converged_reason
5029a09327dSAlp Dener       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
503aad13602SShrirang Abhyankar 
5048e58fa1dSresundermann    test:
5058e58fa1dSresundermann       suffix: 3
5068e58fa1dSresundermann       args: -tao_converged_reason -no_eq
5079a09327dSAlp Dener       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
5088e58fa1dSresundermann 
5098e58fa1dSresundermann    test:
5108e58fa1dSresundermann       suffix: 4
5119a09327dSAlp Dener       args: -tao_converged_reason -tao_almm_type classic
5129a09327dSAlp Dener       requires: !single
5139a09327dSAlp Dener       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
5148e58fa1dSresundermann 
515661095bbSAlp Dener    test:
516661095bbSAlp Dener       suffix: 5
517*c8b0c2b5SHansol Suh       args: -tao_converged_reason -tao_almm_type classic -no_eq
518*c8b0c2b5SHansol Suh       requires: !single
5199a09327dSAlp Dener       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
520661095bbSAlp Dener 
521661095bbSAlp Dener    test:
522661095bbSAlp Dener       suffix: 6
5239a09327dSAlp Dener       args: -tao_converged_reason -tao_almm_subsolver_tao_type bqnktr
5249a09327dSAlp Dener       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
525661095bbSAlp Dener 
526661095bbSAlp Dener    test:
527661095bbSAlp Dener       suffix: 7
5289a09327dSAlp Dener       args: -tao_converged_reason -tao_almm_subsolver_tao_type bncg
5299a09327dSAlp Dener       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
530661095bbSAlp Dener 
531661095bbSAlp Dener    test:
532661095bbSAlp Dener       suffix: 8
533661095bbSAlp Dener       nsize: 2
5349a09327dSAlp Dener       args: -tao_converged_reason
5359a09327dSAlp Dener       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
536661095bbSAlp Dener 
537661095bbSAlp Dener    test:
538661095bbSAlp Dener       suffix: 9
539661095bbSAlp Dener       nsize: 2
5409a09327dSAlp Dener       args: -tao_converged_reason -vec_type cuda -mat_type aijcusparse
5419a09327dSAlp Dener       requires: cuda
5429a09327dSAlp Dener       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
543661095bbSAlp Dener 
544*c8b0c2b5SHansol Suh    test:
545*c8b0c2b5SHansol Suh       suffix: 10
546*c8b0c2b5SHansol Suh       args: -tao_converged_reason -tao_almm_type classic -no_eq -no_bound
547*c8b0c2b5SHansol Suh       requires: !single
548*c8b0c2b5SHansol Suh       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
549*c8b0c2b5SHansol Suh 
550*c8b0c2b5SHansol Suh    test:
551*c8b0c2b5SHansol Suh       suffix: 11
552*c8b0c2b5SHansol Suh       args: -tao_converged_reason -tao_almm_type classic -no_bound
553*c8b0c2b5SHansol Suh       requires: !single
554*c8b0c2b5SHansol Suh       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
555*c8b0c2b5SHansol Suh 
556*c8b0c2b5SHansol Suh    test:
557*c8b0c2b5SHansol Suh       suffix: 12
558*c8b0c2b5SHansol Suh       args: -tao_converged_reason -tao_almm_type phr -no_eq -no_bound
559*c8b0c2b5SHansol Suh       requires: !single
560*c8b0c2b5SHansol Suh       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
561*c8b0c2b5SHansol Suh 
562*c8b0c2b5SHansol Suh    test:
563*c8b0c2b5SHansol Suh       suffix: 13
564*c8b0c2b5SHansol Suh       args: -tao_converged_reason -tao_almm_type phr -no_bound
565*c8b0c2b5SHansol Suh       requires: !single
566*c8b0c2b5SHansol Suh       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"
567*c8b0c2b5SHansol Suh 
568aad13602SShrirang Abhyankar TEST*/
569