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 20aad13602SShrirang Abhyankar static char help[]= "Solves constrained optimiztion 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*f560b561SHong Zhang -init_view : view initial object setup\n\ 25aad13602SShrirang Abhyankar -snes_fd : snes with finite difference Jacobian (needed for pdipm)\n\ 26628da978Sresundermann -snes_compare_explicit : compare user Jacobian with finite difference Jacobian \n\ 27aad13602SShrirang Abhyankar -tao_cmonitor : convergence monitor with constraint norm \n\ 28a5b23f4aSJose E. Roman -tao_view_solution : view exact solution at each iteration\n\ 29*f560b561SHong 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"; 30aad13602SShrirang Abhyankar 31aad13602SShrirang Abhyankar /* 32628da978Sresundermann User-defined application context - contains data needed by the application 33aad13602SShrirang Abhyankar */ 34aad13602SShrirang Abhyankar typedef struct { 35aad13602SShrirang Abhyankar PetscInt n; /* Global length of x */ 36aad13602SShrirang Abhyankar PetscInt ne; /* Global number of equality constraints */ 37aad13602SShrirang Abhyankar PetscInt ni; /* Global number of inequality constraints */ 38*f560b561SHong Zhang PetscBool noeqflag, initview; 39aad13602SShrirang Abhyankar Vec x,xl,xu; 40aad13602SShrirang Abhyankar Vec ce,ci,bl,bu,Xseq; 41aad13602SShrirang Abhyankar Mat Ae,Ai,H; 42aad13602SShrirang Abhyankar VecScatter scat; 43aad13602SShrirang Abhyankar } AppCtx; 44aad13602SShrirang Abhyankar 45aad13602SShrirang Abhyankar /* -------- User-defined Routines --------- */ 46aad13602SShrirang Abhyankar PetscErrorCode InitializeProblem(AppCtx *); 47aad13602SShrirang Abhyankar PetscErrorCode DestroyProblem(AppCtx *); 48aad13602SShrirang Abhyankar PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void *); 49aad13602SShrirang Abhyankar PetscErrorCode FormHessian(Tao,Vec,Mat,Mat, void*); 50aad13602SShrirang Abhyankar PetscErrorCode FormInequalityConstraints(Tao,Vec,Vec,void*); 51aad13602SShrirang Abhyankar PetscErrorCode FormEqualityConstraints(Tao,Vec,Vec,void*); 52aad13602SShrirang Abhyankar PetscErrorCode FormInequalityJacobian(Tao,Vec,Mat,Mat, void*); 53aad13602SShrirang Abhyankar PetscErrorCode FormEqualityJacobian(Tao,Vec,Mat,Mat, void*); 54aad13602SShrirang Abhyankar 55aad13602SShrirang Abhyankar PetscErrorCode main(int argc,char **argv) 56aad13602SShrirang Abhyankar { 57aad13602SShrirang Abhyankar PetscErrorCode ierr; /* used to check for functions returning nonzeros */ 58aad13602SShrirang Abhyankar Tao tao; 59aad13602SShrirang Abhyankar KSP ksp; 60aad13602SShrirang Abhyankar PC pc; 61aad13602SShrirang Abhyankar AppCtx user; /* application context */ 62*f560b561SHong Zhang Vec x,G,CI,CE; 63*f560b561SHong Zhang PetscMPIInt size; 64661095bbSAlp Dener TaoType type; 65*f560b561SHong Zhang PetscReal f; 66661095bbSAlp Dener PetscBool pdipm; 67aad13602SShrirang Abhyankar 68aad13602SShrirang Abhyankar ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 69*f560b561SHong Zhang ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 70*f560b561SHong Zhang if (size>2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_WRONG_MPI_SIZE,"More than 2 processors detected. Example written to use max of 2 processors.\n"); 71aad13602SShrirang Abhyankar 72aad13602SShrirang Abhyankar ierr = PetscPrintf(PETSC_COMM_WORLD,"---- Constrained Problem -----\n");CHKERRQ(ierr); 73aad13602SShrirang Abhyankar ierr = InitializeProblem(&user);CHKERRQ(ierr); /* sets up problem, function below */ 74aad13602SShrirang Abhyankar ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr); 75aad13602SShrirang Abhyankar ierr = TaoSetType(tao,TAOPDIPM);CHKERRQ(ierr); 76aad13602SShrirang Abhyankar ierr = TaoSetInitialVector(tao,user.x);CHKERRQ(ierr); /* gets solution vector from problem */ 77aad13602SShrirang Abhyankar ierr = TaoSetVariableBounds(tao,user.xl,user.xu);CHKERRQ(ierr); /* sets lower upper bounds from given solution */ 78aad13602SShrirang Abhyankar ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void*)&user);CHKERRQ(ierr); 79aad13602SShrirang Abhyankar 808e58fa1dSresundermann if (!user.noeqflag) { 81aad13602SShrirang Abhyankar ierr = TaoSetEqualityConstraintsRoutine(tao,user.ce,FormEqualityConstraints,(void*)&user);CHKERRQ(ierr); 828e58fa1dSresundermann } 83aad13602SShrirang Abhyankar ierr = TaoSetInequalityConstraintsRoutine(tao,user.ci,FormInequalityConstraints,(void*)&user);CHKERRQ(ierr); 848e58fa1dSresundermann if (!user.noeqflag) { 85aad13602SShrirang Abhyankar ierr = TaoSetJacobianEqualityRoutine(tao,user.Ae,user.Ae,FormEqualityJacobian,(void*)&user);CHKERRQ(ierr); /* equality jacobian */ 868e58fa1dSresundermann } 87aad13602SShrirang Abhyankar ierr = TaoSetJacobianInequalityRoutine(tao,user.Ai,user.Ai,FormInequalityJacobian,(void*)&user);CHKERRQ(ierr); /* inequality jacobian */ 88aad13602SShrirang Abhyankar ierr = TaoSetTolerances(tao,1.e-6,1.e-6,1.e-6);CHKERRQ(ierr); 89aad13602SShrirang Abhyankar ierr = TaoSetConstraintTolerances(tao,1.e-6,1.e-6);CHKERRQ(ierr); 90aad13602SShrirang Abhyankar 91aad13602SShrirang Abhyankar ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr); 92aad13602SShrirang Abhyankar ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 9312d688e0SRylee Sundermann ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr); 94aad13602SShrirang Abhyankar /* 95aad13602SShrirang Abhyankar This algorithm produces matrices with zeros along the diagonal therefore we use 9612d688e0SRylee Sundermann MUMPS which provides solver for indefinite matrices 97aad13602SShrirang Abhyankar */ 98*f560b561SHong Zhang #if defined(PETSC_HAVE_MUMPS) 9912d688e0SRylee Sundermann ierr = PCFactorSetMatSolverType(pc,MATSOLVERMUMPS);CHKERRQ(ierr); /* requires mumps to solve pdipm */ 100*f560b561SHong Zhang #else 101*f560b561SHong Zhang if (size > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Requires an external package that supports parallel PCCHOLESKY, e.g., MUMPS."); 102*f560b561SHong Zhang #endif 103aad13602SShrirang Abhyankar ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 104aad13602SShrirang Abhyankar ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 105aad13602SShrirang Abhyankar 106aad13602SShrirang Abhyankar ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 107661095bbSAlp Dener ierr = TaoGetType(tao,&type);CHKERRQ(ierr); 108661095bbSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)tao,TAOPDIPM,&pdipm);CHKERRQ(ierr); 109661095bbSAlp Dener if (pdipm) { 110661095bbSAlp Dener ierr = TaoSetHessianRoutine(tao,user.H,user.H,FormHessian,(void*)&user);CHKERRQ(ierr); 111*f560b561SHong Zhang if (user.initview) { 112*f560b561SHong Zhang ierr = TaoSetUp(tao);CHKERRQ(ierr); 113*f560b561SHong Zhang ierr = VecDuplicate(user.x, &G);CHKERRQ(ierr); 114*f560b561SHong Zhang ierr = FormFunctionGradient(tao, user.x, &f, G, (void*)&user);CHKERRQ(ierr); 115*f560b561SHong Zhang ierr = FormHessian(tao, user.x, user.H, user.H, (void*)&user);CHKERRQ(ierr); 116*f560b561SHong Zhang ierr = PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 117*f560b561SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\nInitial point X:\n",f);CHKERRQ(ierr); 118*f560b561SHong Zhang ierr = VecView(user.x, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 119*f560b561SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\nInitial objective f(x) = %g\n",f);CHKERRQ(ierr); 120*f560b561SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\nInitial gradient and Hessian:\n",f);CHKERRQ(ierr); 121*f560b561SHong Zhang ierr = VecView(G, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 122*f560b561SHong Zhang ierr = MatView(user.H, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 123*f560b561SHong Zhang ierr = VecDestroy(&G);CHKERRQ(ierr); 124*f560b561SHong Zhang ierr = FormInequalityJacobian(tao, user.x, user.Ai, user.Ai, (void*)&user);CHKERRQ(ierr); 125*f560b561SHong Zhang ierr = MatCreateVecs(user.Ai, NULL, &CI);CHKERRQ(ierr); 126*f560b561SHong Zhang ierr = FormInequalityConstraints(tao, user.x, CI, (void*)&user);CHKERRQ(ierr); 127*f560b561SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\nInitial inequality constraints and Jacobian:\n",f);CHKERRQ(ierr); 128*f560b561SHong Zhang ierr = VecView(CI, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 129*f560b561SHong Zhang ierr = MatView(user.Ai, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 130*f560b561SHong Zhang ierr = VecDestroy(&CI);CHKERRQ(ierr); 131*f560b561SHong Zhang if (!user.noeqflag) { 132*f560b561SHong Zhang ierr = FormEqualityJacobian(tao, user.x, user.Ae, user.Ae, (void*)&user);CHKERRQ(ierr); 133*f560b561SHong Zhang ierr = MatCreateVecs(user.Ae, NULL, &CE);CHKERRQ(ierr); 134*f560b561SHong Zhang ierr = FormEqualityConstraints(tao, user.x, CE, (void*)&user);CHKERRQ(ierr); 135*f560b561SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\nInitial equality constraints and Jacobian:\n",f);CHKERRQ(ierr); 136*f560b561SHong Zhang ierr = VecView(CE, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 137*f560b561SHong Zhang ierr = MatView(user.Ae, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 138*f560b561SHong Zhang ierr = VecDestroy(&CE);CHKERRQ(ierr); 139*f560b561SHong Zhang } 140*f560b561SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD, "\n");CHKERRQ(ierr); 141*f560b561SHong Zhang ierr = PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 142*f560b561SHong Zhang } 143661095bbSAlp Dener } 144aad13602SShrirang Abhyankar 145aad13602SShrirang Abhyankar ierr = TaoSolve(tao);CHKERRQ(ierr); 146aad13602SShrirang Abhyankar ierr = TaoGetSolutionVector(tao,&x);CHKERRQ(ierr); 147aad13602SShrirang Abhyankar ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 148aad13602SShrirang Abhyankar 149aad13602SShrirang Abhyankar /* Free objects */ 150aad13602SShrirang Abhyankar ierr = DestroyProblem(&user);CHKERRQ(ierr); 151aad13602SShrirang Abhyankar ierr = TaoDestroy(&tao);CHKERRQ(ierr); 152aad13602SShrirang Abhyankar ierr = PetscFinalize(); 153aad13602SShrirang Abhyankar return ierr; 154aad13602SShrirang Abhyankar } 155aad13602SShrirang Abhyankar 156aad13602SShrirang Abhyankar PetscErrorCode InitializeProblem(AppCtx *user) 157aad13602SShrirang Abhyankar { 158aad13602SShrirang Abhyankar PetscErrorCode ierr; 159aad13602SShrirang Abhyankar PetscMPIInt size; 160aad13602SShrirang Abhyankar PetscMPIInt rank; 161aad13602SShrirang Abhyankar PetscInt nloc,neloc,niloc; 162aad13602SShrirang Abhyankar 163aad13602SShrirang Abhyankar PetscFunctionBegin; 164ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 165ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 1668e58fa1dSresundermann user->noeqflag = PETSC_FALSE; 1678e58fa1dSresundermann ierr = PetscOptionsGetBool(NULL,NULL,"-no_eq",&user->noeqflag,NULL);CHKERRQ(ierr); 168*f560b561SHong Zhang user->initview = PETSC_FALSE; 169*f560b561SHong Zhang ierr = PetscOptionsGetBool(NULL,NULL,"-init_view",&user->initview,NULL);CHKERRQ(ierr); 170628da978Sresundermann 1718e58fa1dSresundermann if (!user->noeqflag) { 172*f560b561SHong Zhang /* Tell user the correct solution, not an error checking */ 1738e58fa1dSresundermann ierr = PetscPrintf(PETSC_COMM_WORLD,"Solution should be f(1,1)=-2\n");CHKERRQ(ierr); 1748e58fa1dSresundermann } 175aad13602SShrirang Abhyankar 1762d4ee042Sprj- /* create vector x and set initial values */ 177aad13602SShrirang Abhyankar user->n = 2; /* global length */ 178*f560b561SHong Zhang nloc = (size==1)?user->n:1; 179661095bbSAlp Dener ierr = VecCreate(PETSC_COMM_WORLD,&user->x);CHKERRQ(ierr); 180661095bbSAlp Dener ierr = VecSetSizes(user->x,nloc,user->n);CHKERRQ(ierr); 181661095bbSAlp Dener ierr = VecSetFromOptions(user->x);CHKERRQ(ierr); 182aad13602SShrirang Abhyankar ierr = VecSet(user->x,0);CHKERRQ(ierr); 183aad13602SShrirang Abhyankar 184aad13602SShrirang Abhyankar /* create and set lower and upper bound vectors */ 185aad13602SShrirang Abhyankar ierr = VecDuplicate(user->x,&user->xl);CHKERRQ(ierr); 186aad13602SShrirang Abhyankar ierr = VecDuplicate(user->x,&user->xu);CHKERRQ(ierr); 187aad13602SShrirang Abhyankar ierr = VecSet(user->xl,-1.0);CHKERRQ(ierr); 188aad13602SShrirang Abhyankar ierr = VecSet(user->xu,2.0);CHKERRQ(ierr); 189aad13602SShrirang Abhyankar 190aad13602SShrirang Abhyankar /* create scater to zero */ 191aad13602SShrirang Abhyankar ierr = VecScatterCreateToZero(user->x,&user->scat,&user->Xseq);CHKERRQ(ierr); 192aad13602SShrirang Abhyankar 193aad13602SShrirang Abhyankar user->ne = 1; 194628da978Sresundermann user->ni = 2; 195aad13602SShrirang Abhyankar neloc = (rank==0)?user->ne:0; 196*f560b561SHong Zhang niloc = (size==1)?user->ni:1; 197aad13602SShrirang Abhyankar 1988e58fa1dSresundermann if (!user->noeqflag) { 199661095bbSAlp Dener ierr = VecCreate(PETSC_COMM_WORLD,&user->ce);CHKERRQ(ierr); /* a 1x1 vec for equality constraints */ 200661095bbSAlp Dener ierr = VecSetSizes(user->ce,neloc,user->ne);CHKERRQ(ierr); 201661095bbSAlp Dener ierr = VecSetFromOptions(user->ce);CHKERRQ(ierr); 202661095bbSAlp Dener ierr = VecSetUp(user->ce);CHKERRQ(ierr); 2038e58fa1dSresundermann } 204628da978Sresundermann 205661095bbSAlp Dener ierr = VecCreate(PETSC_COMM_WORLD,&user->ci);CHKERRQ(ierr); /* a 2x1 vec for inequality constraints */ 206661095bbSAlp Dener ierr = VecSetSizes(user->ci,niloc,user->ni);CHKERRQ(ierr); 207661095bbSAlp Dener ierr = VecSetFromOptions(user->ci);CHKERRQ(ierr); 208661095bbSAlp Dener ierr = VecSetUp(user->ci);CHKERRQ(ierr); 209aad13602SShrirang Abhyankar 210a5b23f4aSJose E. Roman /* nexn & nixn matricies for equally and inequalty constraints */ 2118e58fa1dSresundermann if (!user->noeqflag) { 212aad13602SShrirang Abhyankar ierr = MatCreate(PETSC_COMM_WORLD,&user->Ae);CHKERRQ(ierr); 2138e58fa1dSresundermann ierr = MatSetSizes(user->Ae,neloc,nloc,user->ne,user->n);CHKERRQ(ierr); 2148e58fa1dSresundermann ierr = MatSetFromOptions(user->Ae);CHKERRQ(ierr); 215661095bbSAlp Dener ierr = MatSetUp(user->Ae);CHKERRQ(ierr); 2168e58fa1dSresundermann } 2178e58fa1dSresundermann 218aad13602SShrirang Abhyankar ierr = MatCreate(PETSC_COMM_WORLD,&user->Ai);CHKERRQ(ierr); 219aad13602SShrirang Abhyankar ierr = MatSetSizes(user->Ai,niloc,nloc,user->ni,user->n);CHKERRQ(ierr); 220aad13602SShrirang Abhyankar ierr = MatSetFromOptions(user->Ai);CHKERRQ(ierr); 221661095bbSAlp Dener ierr = MatSetUp(user->Ai);CHKERRQ(ierr); 222628da978Sresundermann 223628da978Sresundermann ierr = MatCreate(PETSC_COMM_WORLD,&user->H);CHKERRQ(ierr); 224628da978Sresundermann ierr = MatSetSizes(user->H,nloc,nloc,user->n,user->n);CHKERRQ(ierr); 225628da978Sresundermann ierr = MatSetFromOptions(user->H);CHKERRQ(ierr); 226661095bbSAlp Dener ierr = MatSetUp(user->H);CHKERRQ(ierr); 227aad13602SShrirang Abhyankar PetscFunctionReturn(0); 228aad13602SShrirang Abhyankar } 229aad13602SShrirang Abhyankar 230aad13602SShrirang Abhyankar PetscErrorCode DestroyProblem(AppCtx *user) 231aad13602SShrirang Abhyankar { 232aad13602SShrirang Abhyankar PetscErrorCode ierr; 233aad13602SShrirang Abhyankar 234aad13602SShrirang Abhyankar PetscFunctionBegin; 2358e58fa1dSresundermann if (!user->noeqflag) { 236aad13602SShrirang Abhyankar ierr = MatDestroy(&user->Ae);CHKERRQ(ierr); 2378e58fa1dSresundermann } 238aad13602SShrirang Abhyankar ierr = MatDestroy(&user->Ai);CHKERRQ(ierr); 239aad13602SShrirang Abhyankar ierr = MatDestroy(&user->H);CHKERRQ(ierr); 240aad13602SShrirang Abhyankar 241aad13602SShrirang Abhyankar ierr = VecDestroy(&user->x);CHKERRQ(ierr); 2428e58fa1dSresundermann if (!user->noeqflag) { 243aad13602SShrirang Abhyankar ierr = VecDestroy(&user->ce);CHKERRQ(ierr); 2448e58fa1dSresundermann } 245aad13602SShrirang Abhyankar ierr = VecDestroy(&user->ci);CHKERRQ(ierr); 246aad13602SShrirang Abhyankar ierr = VecDestroy(&user->xl);CHKERRQ(ierr); 247aad13602SShrirang Abhyankar ierr = VecDestroy(&user->xu);CHKERRQ(ierr); 248aad13602SShrirang Abhyankar ierr = VecDestroy(&user->Xseq);CHKERRQ(ierr); 249aad13602SShrirang Abhyankar ierr = VecScatterDestroy(&user->scat);CHKERRQ(ierr); 250aad13602SShrirang Abhyankar PetscFunctionReturn(0); 251aad13602SShrirang Abhyankar } 252aad13602SShrirang Abhyankar 253628da978Sresundermann /* Evaluate 254628da978Sresundermann f(x) = (x0 - 2)^2 + (x1 - 2)^2 - 2*(x0 + x1) 255628da978Sresundermann G = grad f = [2*(x0 - 2) - 2; 256628da978Sresundermann 2*(x1 - 2) - 2] 257aad13602SShrirang Abhyankar */ 258aad13602SShrirang Abhyankar PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ctx) 259aad13602SShrirang Abhyankar { 260aad13602SShrirang Abhyankar PetscScalar g; 261aad13602SShrirang Abhyankar const PetscScalar *x; 262aad13602SShrirang Abhyankar MPI_Comm comm; 263aad13602SShrirang Abhyankar PetscMPIInt rank; 264aad13602SShrirang Abhyankar PetscErrorCode ierr; 265aad13602SShrirang Abhyankar PetscReal fin; 266aad13602SShrirang Abhyankar AppCtx *user=(AppCtx*)ctx; 267aad13602SShrirang Abhyankar Vec Xseq=user->Xseq; 268aad13602SShrirang Abhyankar VecScatter scat=user->scat; 269aad13602SShrirang Abhyankar 270aad13602SShrirang Abhyankar PetscFunctionBegin; 271aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 272ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 273aad13602SShrirang Abhyankar 274aad13602SShrirang Abhyankar ierr = VecScatterBegin(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 275aad13602SShrirang Abhyankar ierr = VecScatterEnd(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 276aad13602SShrirang Abhyankar 277aad13602SShrirang Abhyankar fin = 0.0; 278aad13602SShrirang Abhyankar if (!rank) { 279aad13602SShrirang Abhyankar ierr = VecGetArrayRead(Xseq,&x);CHKERRQ(ierr); 280aad13602SShrirang Abhyankar fin = (x[0]-2.0)*(x[0]-2.0) + (x[1]-2.0)*(x[1]-2.0) - 2.0*(x[0]+x[1]); 281aad13602SShrirang Abhyankar g = 2.0*(x[0]-2.0) - 2.0; 282aad13602SShrirang Abhyankar ierr = VecSetValue(G,0,g,INSERT_VALUES);CHKERRQ(ierr); 283aad13602SShrirang Abhyankar g = 2.0*(x[1]-2.0) - 2.0; 284aad13602SShrirang Abhyankar ierr = VecSetValue(G,1,g,INSERT_VALUES);CHKERRQ(ierr); 285aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(Xseq,&x);CHKERRQ(ierr); 286aad13602SShrirang Abhyankar } 287ffc4695bSBarry Smith ierr = MPI_Allreduce(&fin,f,1,MPIU_REAL,MPIU_SUM,comm);CHKERRMPI(ierr); 288aad13602SShrirang Abhyankar ierr = VecAssemblyBegin(G);CHKERRQ(ierr); 289aad13602SShrirang Abhyankar ierr = VecAssemblyEnd(G);CHKERRQ(ierr); 290aad13602SShrirang Abhyankar PetscFunctionReturn(0); 291aad13602SShrirang Abhyankar } 292aad13602SShrirang Abhyankar 293628da978Sresundermann /* Evaluate 294628da978Sresundermann H = fxx + grad (grad g^T*DI) - grad (grad h^T*DE)] 295628da978Sresundermann = [ 2*(1+de[0]-di[0]+di[1]), 0; 296628da978Sresundermann 0, 2] 297628da978Sresundermann */ 298aad13602SShrirang Abhyankar PetscErrorCode FormHessian(Tao tao, Vec x,Mat H, Mat Hpre, void *ctx) 299aad13602SShrirang Abhyankar { 300aad13602SShrirang Abhyankar AppCtx *user=(AppCtx*)ctx; 301aad13602SShrirang Abhyankar Vec DE,DI; 302aad13602SShrirang Abhyankar const PetscScalar *de,*di; 303aad13602SShrirang Abhyankar PetscInt zero=0,one=1; 304aad13602SShrirang Abhyankar PetscScalar two=2.0; 305aad13602SShrirang Abhyankar PetscScalar val=0.0; 306*f560b561SHong Zhang Vec Deseq,Diseq; 307*f560b561SHong Zhang VecScatter Descat,Discat; 308aad13602SShrirang Abhyankar PetscMPIInt rank; 309aad13602SShrirang Abhyankar MPI_Comm comm; 310aad13602SShrirang Abhyankar PetscErrorCode ierr; 311aad13602SShrirang Abhyankar 312aad13602SShrirang Abhyankar PetscFunctionBegin; 313aad13602SShrirang Abhyankar ierr = TaoGetDualVariables(tao,&DE,&DI);CHKERRQ(ierr); 314aad13602SShrirang Abhyankar 315aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 316ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 317aad13602SShrirang Abhyankar 3188e58fa1dSresundermann if (!user->noeqflag) { 319aad13602SShrirang Abhyankar ierr = VecScatterCreateToZero(DE,&Descat,&Deseq);CHKERRQ(ierr); 320aad13602SShrirang Abhyankar ierr = VecScatterBegin(Descat,DE,Deseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 321aad13602SShrirang Abhyankar ierr = VecScatterEnd(Descat,DE,Deseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3228e58fa1dSresundermann } 323*f560b561SHong Zhang ierr = VecScatterCreateToZero(DI,&Discat,&Diseq);CHKERRQ(ierr); 3248e58fa1dSresundermann ierr = VecScatterBegin(Discat,DI,Diseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 325aad13602SShrirang Abhyankar ierr = VecScatterEnd(Discat,DI,Diseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 326aad13602SShrirang Abhyankar 327aad13602SShrirang Abhyankar if (!rank) { 3288e58fa1dSresundermann if (!user->noeqflag) { 329aad13602SShrirang Abhyankar ierr = VecGetArrayRead(Deseq,&de);CHKERRQ(ierr); /* places equality constraint dual into array */ 3308e58fa1dSresundermann } 331aad13602SShrirang Abhyankar ierr = VecGetArrayRead(Diseq,&di);CHKERRQ(ierr); /* places inequality constraint dual into array */ 332628da978Sresundermann 3338e58fa1dSresundermann if (!user->noeqflag) { 334628da978Sresundermann val = 2.0 * (1 + de[0] - di[0] + di[1]); 335aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(Deseq,&de);CHKERRQ(ierr); 336628da978Sresundermann ierr = VecRestoreArrayRead(Diseq,&di);CHKERRQ(ierr); 3378e58fa1dSresundermann } else { 338628da978Sresundermann val = 2.0 * (1 - di[0] + di[1]); 3398e58fa1dSresundermann } 340aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(Diseq,&di);CHKERRQ(ierr); 341aad13602SShrirang Abhyankar ierr = MatSetValues(H,1,&zero,1,&zero,&val,INSERT_VALUES);CHKERRQ(ierr); 342aad13602SShrirang Abhyankar ierr = MatSetValues(H,1,&one,1,&one,&two,INSERT_VALUES);CHKERRQ(ierr); 343aad13602SShrirang Abhyankar } 344aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 345aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3468e58fa1dSresundermann if (!user->noeqflag) { 347aad13602SShrirang Abhyankar ierr = VecScatterDestroy(&Descat);CHKERRQ(ierr); 348aad13602SShrirang Abhyankar ierr = VecDestroy(&Deseq);CHKERRQ(ierr); 3498e58fa1dSresundermann } 350*f560b561SHong Zhang ierr = VecScatterDestroy(&Discat);CHKERRQ(ierr); 351*f560b561SHong Zhang ierr = VecDestroy(&Diseq);CHKERRQ(ierr); 352aad13602SShrirang Abhyankar PetscFunctionReturn(0); 353aad13602SShrirang Abhyankar } 354aad13602SShrirang Abhyankar 355628da978Sresundermann /* Evaluate 356628da978Sresundermann h = [ x0^2 - x1; 357628da978Sresundermann 1 -(x0^2 - x1)] 358628da978Sresundermann */ 359aad13602SShrirang Abhyankar PetscErrorCode FormInequalityConstraints(Tao tao,Vec X,Vec CI,void *ctx) 360aad13602SShrirang Abhyankar { 361aad13602SShrirang Abhyankar const PetscScalar *x; 362aad13602SShrirang Abhyankar PetscScalar ci; 363aad13602SShrirang Abhyankar PetscErrorCode ierr; 364aad13602SShrirang Abhyankar MPI_Comm comm; 365aad13602SShrirang Abhyankar PetscMPIInt rank; 366aad13602SShrirang Abhyankar AppCtx *user=(AppCtx*)ctx; 367aad13602SShrirang Abhyankar Vec Xseq=user->Xseq; 368aad13602SShrirang Abhyankar VecScatter scat=user->scat; 369aad13602SShrirang Abhyankar 370aad13602SShrirang Abhyankar PetscFunctionBegin; 371aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 372ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 373aad13602SShrirang Abhyankar 374aad13602SShrirang Abhyankar ierr = VecScatterBegin(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 375aad13602SShrirang Abhyankar ierr = VecScatterEnd(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 376aad13602SShrirang Abhyankar 377aad13602SShrirang Abhyankar if (!rank) { 378aad13602SShrirang Abhyankar ierr = VecGetArrayRead(Xseq,&x);CHKERRQ(ierr); 379aad13602SShrirang Abhyankar ci = x[0]*x[0] - x[1]; 380aad13602SShrirang Abhyankar ierr = VecSetValue(CI,0,ci,INSERT_VALUES);CHKERRQ(ierr); 381aad13602SShrirang Abhyankar ci = -x[0]*x[0] + x[1] + 1.0; 382aad13602SShrirang Abhyankar ierr = VecSetValue(CI,1,ci,INSERT_VALUES);CHKERRQ(ierr); 383aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(Xseq,&x);CHKERRQ(ierr); 384aad13602SShrirang Abhyankar } 385aad13602SShrirang Abhyankar ierr = VecAssemblyBegin(CI);CHKERRQ(ierr); 386aad13602SShrirang Abhyankar ierr = VecAssemblyEnd(CI);CHKERRQ(ierr); 387aad13602SShrirang Abhyankar PetscFunctionReturn(0); 388aad13602SShrirang Abhyankar } 389aad13602SShrirang Abhyankar 390628da978Sresundermann /* Evaluate 391628da978Sresundermann g = [ x0^2 + x1 - 2] 392628da978Sresundermann */ 393aad13602SShrirang Abhyankar PetscErrorCode FormEqualityConstraints(Tao tao,Vec X,Vec CE,void *ctx) 394aad13602SShrirang Abhyankar { 395aad13602SShrirang Abhyankar const PetscScalar *x; 396aad13602SShrirang Abhyankar PetscScalar ce; 397aad13602SShrirang Abhyankar PetscErrorCode ierr; 398aad13602SShrirang Abhyankar MPI_Comm comm; 399aad13602SShrirang Abhyankar PetscMPIInt rank; 400aad13602SShrirang Abhyankar AppCtx *user=(AppCtx*)ctx; 401aad13602SShrirang Abhyankar Vec Xseq=user->Xseq; 402aad13602SShrirang Abhyankar VecScatter scat=user->scat; 403aad13602SShrirang Abhyankar 404aad13602SShrirang Abhyankar PetscFunctionBegin; 405aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 406ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 407aad13602SShrirang Abhyankar 408aad13602SShrirang Abhyankar ierr = VecScatterBegin(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 409aad13602SShrirang Abhyankar ierr = VecScatterEnd(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 410aad13602SShrirang Abhyankar 411aad13602SShrirang Abhyankar if (!rank) { 412aad13602SShrirang Abhyankar ierr = VecGetArrayRead(Xseq,&x);CHKERRQ(ierr); 413aad13602SShrirang Abhyankar ce = x[0]*x[0] + x[1] - 2.0; 414aad13602SShrirang Abhyankar ierr = VecSetValue(CE,0,ce,INSERT_VALUES);CHKERRQ(ierr); 415aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(Xseq,&x);CHKERRQ(ierr); 416aad13602SShrirang Abhyankar } 417aad13602SShrirang Abhyankar ierr = VecAssemblyBegin(CE);CHKERRQ(ierr); 418aad13602SShrirang Abhyankar ierr = VecAssemblyEnd(CE);CHKERRQ(ierr); 419aad13602SShrirang Abhyankar PetscFunctionReturn(0); 420aad13602SShrirang Abhyankar } 421aad13602SShrirang Abhyankar 422628da978Sresundermann /* 423628da978Sresundermann grad h = [ 2*x0, -1; 424628da978Sresundermann -2*x0, 1] 425628da978Sresundermann */ 426aad13602SShrirang Abhyankar PetscErrorCode FormInequalityJacobian(Tao tao, Vec X, Mat JI, Mat JIpre, void *ctx) 427aad13602SShrirang Abhyankar { 428aad13602SShrirang Abhyankar AppCtx *user=(AppCtx*)ctx; 429*f560b561SHong Zhang PetscInt zero=0,one=1,cols[2]; 430aad13602SShrirang Abhyankar PetscScalar vals[2]; 431aad13602SShrirang Abhyankar const PetscScalar *x; 432aad13602SShrirang Abhyankar PetscErrorCode ierr; 433aad13602SShrirang Abhyankar Vec Xseq=user->Xseq; 434aad13602SShrirang Abhyankar VecScatter scat=user->scat; 435aad13602SShrirang Abhyankar MPI_Comm comm; 436aad13602SShrirang Abhyankar PetscMPIInt rank; 437aad13602SShrirang Abhyankar 438aad13602SShrirang Abhyankar PetscFunctionBegin; 439aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 440ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 441aad13602SShrirang Abhyankar ierr = VecScatterBegin(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 442aad13602SShrirang Abhyankar ierr = VecScatterEnd(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 443aad13602SShrirang Abhyankar 444aad13602SShrirang Abhyankar ierr = VecGetArrayRead(Xseq,&x);CHKERRQ(ierr); 445*f560b561SHong Zhang if (!rank) { 446aad13602SShrirang Abhyankar cols[0] = 0; cols[1] = 1; 447628da978Sresundermann vals[0] = 2*x[0]; vals[1] = -1.0; 448*f560b561SHong Zhang ierr = MatSetValues(JI,1,&zero,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 449628da978Sresundermann vals[0] = -2*x[0]; vals[1] = 1.0; 450*f560b561SHong Zhang ierr = MatSetValues(JI,1,&one,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 451aad13602SShrirang Abhyankar } 452aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(Xseq,&x);CHKERRQ(ierr); 453aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(JI,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 454aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(JI,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 455aad13602SShrirang Abhyankar PetscFunctionReturn(0); 456aad13602SShrirang Abhyankar } 457aad13602SShrirang Abhyankar 458628da978Sresundermann /* 459628da978Sresundermann grad g = [2*x0 460628da978Sresundermann 1.0 ] 461628da978Sresundermann */ 462aad13602SShrirang Abhyankar PetscErrorCode FormEqualityJacobian(Tao tao,Vec X,Mat JE,Mat JEpre,void *ctx) 463aad13602SShrirang Abhyankar { 464*f560b561SHong Zhang PetscInt zero=0,cols[2]; 465aad13602SShrirang Abhyankar PetscScalar vals[2]; 466aad13602SShrirang Abhyankar const PetscScalar *x; 467aad13602SShrirang Abhyankar PetscMPIInt rank; 468aad13602SShrirang Abhyankar MPI_Comm comm; 469aad13602SShrirang Abhyankar PetscErrorCode ierr; 470aad13602SShrirang Abhyankar 471aad13602SShrirang Abhyankar PetscFunctionBegin; 472aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 473ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 474aad13602SShrirang Abhyankar 475aad13602SShrirang Abhyankar if (!rank) { 476aad13602SShrirang Abhyankar ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 477*f560b561SHong Zhang cols[0] = 0; cols[1] = 1; 478aad13602SShrirang Abhyankar vals[0] = 2*x[0]; vals[1] = 1.0; 479*f560b561SHong Zhang ierr = MatSetValues(JE,1,&zero,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 480aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 481aad13602SShrirang Abhyankar } 482aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(JE,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 483aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(JE,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 484aad13602SShrirang Abhyankar PetscFunctionReturn(0); 485aad13602SShrirang Abhyankar } 486aad13602SShrirang Abhyankar 487aad13602SShrirang Abhyankar /*TEST 488aad13602SShrirang Abhyankar 489aad13602SShrirang Abhyankar build: 490*f560b561SHong Zhang requires: !complex !defined(PETSC_USE_CXX) 491aad13602SShrirang Abhyankar 492aad13602SShrirang Abhyankar test: 49309ee8bb0SRylee Sundermann args: -tao_converged_reason -tao_pdipm_kkt_shift_pd 494*f560b561SHong Zhang requires: mumps 495aad13602SShrirang Abhyankar 496aad13602SShrirang Abhyankar test: 497aad13602SShrirang Abhyankar suffix: 2 498aad13602SShrirang Abhyankar nsize: 2 49909ee8bb0SRylee Sundermann args: -tao_converged_reason -tao_pdipm_kkt_shift_pd 500*f560b561SHong Zhang requires: mumps 501aad13602SShrirang Abhyankar 5028e58fa1dSresundermann test: 5038e58fa1dSresundermann suffix: 3 5048e58fa1dSresundermann args: -tao_converged_reason -no_eq 505*f560b561SHong Zhang requires: mumps 5068e58fa1dSresundermann 5078e58fa1dSresundermann test: 5088e58fa1dSresundermann suffix: 4 5098e58fa1dSresundermann nsize: 2 5108e58fa1dSresundermann args: -tao_converged_reason -no_eq 511*f560b561SHong Zhang requires: mumps 5128e58fa1dSresundermann 513661095bbSAlp Dener test: 514661095bbSAlp Dener suffix: 5 515661095bbSAlp Dener args: -tao_cmonitor -tao_type almm 516*f560b561SHong Zhang requires: mumps 517661095bbSAlp Dener 518661095bbSAlp Dener test: 519661095bbSAlp Dener suffix: 6 520661095bbSAlp Dener args: -tao_cmonitor -tao_type almm -tao_almm_type phr 521*f560b561SHong Zhang requires: mumps 522661095bbSAlp Dener 523661095bbSAlp Dener test: 524661095bbSAlp Dener suffix: 7 525661095bbSAlp Dener nsize: 2 526*f560b561SHong Zhang requires: mumps 527661095bbSAlp Dener args: -tao_cmonitor -tao_type almm 528661095bbSAlp Dener 529661095bbSAlp Dener test: 530661095bbSAlp Dener suffix: 8 531661095bbSAlp Dener nsize: 2 532*f560b561SHong Zhang requires: cuda mumps 533661095bbSAlp Dener args: -tao_cmonitor -tao_type almm -vec_type cuda -mat_type aijcusparse 534661095bbSAlp Dener 535661095bbSAlp Dener test: 536661095bbSAlp Dener suffix: 9 537661095bbSAlp Dener nsize: 2 538661095bbSAlp Dener args: -tao_cmonitor -tao_type almm -no_eq 539*f560b561SHong Zhang requires: mumps 540661095bbSAlp Dener 541661095bbSAlp Dener test: 542661095bbSAlp Dener suffix: 10 543661095bbSAlp Dener nsize: 2 544661095bbSAlp Dener args: -tao_cmonitor -tao_type almm -tao_almm_type phr -no_eq 545*f560b561SHong Zhang requires: mumps 546661095bbSAlp Dener 547aad13602SShrirang Abhyankar TEST*/ 548