1aad13602SShrirang Abhyankar /* Program usage: mpiexec -n 2 ./ex1 [-help] [all TAO options] */ 2aad13602SShrirang Abhyankar 3aad13602SShrirang Abhyankar /* ---------------------------------------------------------------------- 4aad13602SShrirang Abhyankar min f = (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 8aad13602SShrirang Abhyankar ---------------------------------------------------------------------- */ 9aad13602SShrirang Abhyankar 10aad13602SShrirang Abhyankar #include <petsctao.h> 11aad13602SShrirang Abhyankar 12aad13602SShrirang Abhyankar static char help[]= "Solves constrained optimiztion problem using pdipm.\n\ 13aad13602SShrirang Abhyankar Input parameters include:\n\ 14aad13602SShrirang Abhyankar -tao_type pdipm : sets Tao solver\n\ 15*8e58fa1dSresundermann -no_eq : removes the equaility constraints from the problem\n\ 16aad13602SShrirang Abhyankar -snes_fd : snes with finite difference Jacobian (needed for pdipm)\n\ 17aad13602SShrirang Abhyankar -tao_cmonitor : convergence monitor with constraint norm \n\ 18aad13602SShrirang Abhyankar -tao_view_solution : view exact solution at each itteration\n\ 19aad13602SShrirang Abhyankar Note: external package superlu_dist is requried to run either for ipm or pdipm. Additionally This is designed for a maximum of 2 processors, the code will error if size > 2.\n\n"; 20aad13602SShrirang Abhyankar 21aad13602SShrirang Abhyankar /* 22aad13602SShrirang Abhyankar User-defined application context - contains data needed by the 23aad13602SShrirang Abhyankar application-provided call-back routines, FormFunction(), 24aad13602SShrirang Abhyankar FormGradient(), and FormHessian(). 25aad13602SShrirang Abhyankar */ 26aad13602SShrirang Abhyankar 27aad13602SShrirang Abhyankar /* 28aad13602SShrirang Abhyankar x,d in R^n 29aad13602SShrirang Abhyankar f in R 30aad13602SShrirang Abhyankar bin in R^mi 31aad13602SShrirang Abhyankar beq in R^me 32aad13602SShrirang Abhyankar Aeq in R^(me x n) 33aad13602SShrirang Abhyankar Ain in R^(mi x n) 34aad13602SShrirang Abhyankar H in R^(n x n) 35aad13602SShrirang Abhyankar min f=(1/2)*x'*H*x + d'*x 36aad13602SShrirang Abhyankar s.t. Aeq*x == beq 37aad13602SShrirang Abhyankar Ain*x >= bin 38aad13602SShrirang Abhyankar */ 39aad13602SShrirang Abhyankar typedef struct { 40aad13602SShrirang Abhyankar PetscInt n; /* Global length of x */ 41aad13602SShrirang Abhyankar PetscInt ne; /* Global number of equality constraints */ 42aad13602SShrirang Abhyankar PetscInt ni; /* Global number of inequality constraints */ 43*8e58fa1dSresundermann PetscBool noeqflag; 44aad13602SShrirang Abhyankar Vec x,xl,xu; 45aad13602SShrirang Abhyankar Vec ce,ci,bl,bu,Xseq; 46aad13602SShrirang Abhyankar Mat Ae,Ai,H; 47aad13602SShrirang Abhyankar VecScatter scat; 48aad13602SShrirang Abhyankar } AppCtx; 49aad13602SShrirang Abhyankar 50aad13602SShrirang Abhyankar /* -------- User-defined Routines --------- */ 51aad13602SShrirang Abhyankar PetscErrorCode InitializeProblem(AppCtx *); 52aad13602SShrirang Abhyankar PetscErrorCode DestroyProblem(AppCtx *); 53aad13602SShrirang Abhyankar PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void *); 54aad13602SShrirang Abhyankar PetscErrorCode FormHessian(Tao,Vec,Mat,Mat, void*); 55aad13602SShrirang Abhyankar PetscErrorCode FormInequalityConstraints(Tao,Vec,Vec,void*); 56aad13602SShrirang Abhyankar PetscErrorCode FormEqualityConstraints(Tao,Vec,Vec,void*); 57aad13602SShrirang Abhyankar PetscErrorCode FormInequalityJacobian(Tao,Vec,Mat,Mat, void*); 58aad13602SShrirang Abhyankar PetscErrorCode FormEqualityJacobian(Tao,Vec,Mat,Mat, void*); 59aad13602SShrirang Abhyankar 60aad13602SShrirang Abhyankar PetscErrorCode main(int argc,char **argv) 61aad13602SShrirang Abhyankar { 62aad13602SShrirang Abhyankar PetscErrorCode ierr; /* used to check for functions returning nonzeros */ 63aad13602SShrirang Abhyankar Tao tao; 64aad13602SShrirang Abhyankar KSP ksp; 65aad13602SShrirang Abhyankar PC pc; 66aad13602SShrirang Abhyankar AppCtx user; /* application context */ 67aad13602SShrirang Abhyankar Vec x; 68aad13602SShrirang Abhyankar PetscMPIInt rank; 69aad13602SShrirang Abhyankar 70aad13602SShrirang Abhyankar ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 71aad13602SShrirang Abhyankar ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 72aad13602SShrirang Abhyankar if (rank>1){ 73aad13602SShrirang Abhyankar SETERRQ(PETSC_COMM_SELF,1,"More than 2 processors detected. Example written to use max of 2 processors.\n"); 74aad13602SShrirang Abhyankar } 75aad13602SShrirang Abhyankar 76aad13602SShrirang Abhyankar ierr = PetscPrintf(PETSC_COMM_WORLD,"---- Constrained Problem -----\n");CHKERRQ(ierr); 77aad13602SShrirang Abhyankar ierr = InitializeProblem(&user);CHKERRQ(ierr); /* sets up problem, function below */ 78aad13602SShrirang Abhyankar ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr); 79aad13602SShrirang Abhyankar ierr = TaoSetType(tao,TAOPDIPM);CHKERRQ(ierr); 80aad13602SShrirang Abhyankar ierr = TaoSetInitialVector(tao,user.x);CHKERRQ(ierr); /* gets solution vector from problem */ 81aad13602SShrirang Abhyankar ierr = TaoSetVariableBounds(tao,user.xl,user.xu);CHKERRQ(ierr); /* sets lower upper bounds from given solution */ 82aad13602SShrirang Abhyankar ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void*)&user);CHKERRQ(ierr); 83aad13602SShrirang Abhyankar 84*8e58fa1dSresundermann if (!user.noeqflag){ 85aad13602SShrirang Abhyankar ierr = TaoSetEqualityConstraintsRoutine(tao,user.ce,FormEqualityConstraints,(void*)&user);CHKERRQ(ierr); 86*8e58fa1dSresundermann } 87aad13602SShrirang Abhyankar ierr = TaoSetInequalityConstraintsRoutine(tao,user.ci,FormInequalityConstraints,(void*)&user);CHKERRQ(ierr); 88aad13602SShrirang Abhyankar 89*8e58fa1dSresundermann if (!user.noeqflag){ 90aad13602SShrirang Abhyankar ierr = TaoSetJacobianEqualityRoutine(tao,user.Ae,user.Ae,FormEqualityJacobian,(void*)&user);CHKERRQ(ierr); /* equality jacobian */ 91*8e58fa1dSresundermann } 92aad13602SShrirang Abhyankar ierr = TaoSetJacobianInequalityRoutine(tao,user.Ai,user.Ai,FormInequalityJacobian,(void*)&user);CHKERRQ(ierr); /* inequality jacobian */ 93aad13602SShrirang Abhyankar ierr = TaoSetHessianRoutine(tao,user.H,user.H,FormHessian,(void*)&user);CHKERRQ(ierr); 94aad13602SShrirang Abhyankar ierr = TaoSetTolerances(tao,1.e-6,1.e-6,1.e-6);CHKERRQ(ierr); 95aad13602SShrirang Abhyankar ierr = TaoSetConstraintTolerances(tao,1.e-6,1.e-6);CHKERRQ(ierr); 96aad13602SShrirang Abhyankar 97aad13602SShrirang Abhyankar ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr); 98aad13602SShrirang Abhyankar ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 99aad13602SShrirang Abhyankar ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 100aad13602SShrirang Abhyankar /* 101aad13602SShrirang Abhyankar This algorithm produces matrices with zeros along the diagonal therefore we use 102aad13602SShrirang Abhyankar SuperLU_DIST which provides shift to the diagonal 103aad13602SShrirang Abhyankar */ 104aad13602SShrirang Abhyankar ierr = PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr); /* requires superlu_dist to solve pdipm */ 105aad13602SShrirang Abhyankar ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 106aad13602SShrirang Abhyankar ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 107aad13602SShrirang Abhyankar 108aad13602SShrirang Abhyankar ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 109aad13602SShrirang Abhyankar 110aad13602SShrirang Abhyankar ierr = TaoSolve(tao);CHKERRQ(ierr); 111aad13602SShrirang Abhyankar ierr = TaoGetSolutionVector(tao,&x);CHKERRQ(ierr); 112aad13602SShrirang Abhyankar ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 113aad13602SShrirang Abhyankar 114aad13602SShrirang Abhyankar /* Free objects */ 115aad13602SShrirang Abhyankar ierr = DestroyProblem(&user);CHKERRQ(ierr); 116aad13602SShrirang Abhyankar ierr = TaoDestroy(&tao);CHKERRQ(ierr); 117aad13602SShrirang Abhyankar ierr = PetscFinalize(); 118aad13602SShrirang Abhyankar return ierr; 119aad13602SShrirang Abhyankar } 120aad13602SShrirang Abhyankar 121aad13602SShrirang Abhyankar PetscErrorCode InitializeProblem(AppCtx *user) 122aad13602SShrirang Abhyankar { 123aad13602SShrirang Abhyankar PetscErrorCode ierr; 124aad13602SShrirang Abhyankar PetscMPIInt size; 125aad13602SShrirang Abhyankar PetscMPIInt rank; 126aad13602SShrirang Abhyankar PetscInt nloc,neloc,niloc; 127aad13602SShrirang Abhyankar 128aad13602SShrirang Abhyankar PetscFunctionBegin; 129aad13602SShrirang Abhyankar ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 130aad13602SShrirang Abhyankar ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 131*8e58fa1dSresundermann user->noeqflag = PETSC_FALSE; 132*8e58fa1dSresundermann ierr = PetscOptionsGetBool(NULL,NULL,"-no_eq",&user->noeqflag,NULL);CHKERRQ(ierr); 133*8e58fa1dSresundermann if (!user->noeqflag) { 134*8e58fa1dSresundermann ierr = PetscPrintf(PETSC_COMM_WORLD,"Solution should be f(1,1)=-2\n");CHKERRQ(ierr); 135*8e58fa1dSresundermann } 136aad13602SShrirang Abhyankar 1372d4ee042Sprj- /* create vector x and set initial values */ 138aad13602SShrirang Abhyankar user->n = 2; /* global length */ 139aad13602SShrirang Abhyankar nloc = (rank==0)?user->n:0; 140aad13602SShrirang Abhyankar ierr = VecCreateMPI(PETSC_COMM_WORLD,nloc,user->n,&user->x);CHKERRQ(ierr); 141aad13602SShrirang Abhyankar ierr = VecSet(user->x,0);CHKERRQ(ierr); 142aad13602SShrirang Abhyankar 143aad13602SShrirang Abhyankar /* create and set lower and upper bound vectors */ 144aad13602SShrirang Abhyankar ierr = VecDuplicate(user->x,&user->xl);CHKERRQ(ierr); 145aad13602SShrirang Abhyankar ierr = VecDuplicate(user->x,&user->xu);CHKERRQ(ierr); 146aad13602SShrirang Abhyankar ierr = VecSet(user->xl,-1.0);CHKERRQ(ierr); 147aad13602SShrirang Abhyankar ierr = VecSet(user->xu,2.0);CHKERRQ(ierr); 148aad13602SShrirang Abhyankar 149aad13602SShrirang Abhyankar /* create scater to zero */ 150aad13602SShrirang Abhyankar ierr = VecScatterCreateToZero(user->x,&user->scat,&user->Xseq);CHKERRQ(ierr); 151aad13602SShrirang Abhyankar 152aad13602SShrirang Abhyankar user->ne = 1; 153aad13602SShrirang Abhyankar neloc = (rank==0)?user->ne:0; 154aad13602SShrirang Abhyankar 155*8e58fa1dSresundermann if (!user->noeqflag){ 156*8e58fa1dSresundermann ierr = VecCreateMPI(PETSC_COMM_WORLD,neloc,user->ne,&user->ce);CHKERRQ(ierr); /* a 1x1 vec for equality constraints */ 157*8e58fa1dSresundermann } 158aad13602SShrirang Abhyankar user->ni = 2; 159aad13602SShrirang Abhyankar niloc = (rank==0)?user->ni:0; 160aad13602SShrirang Abhyankar ierr = VecCreateMPI(PETSC_COMM_WORLD,niloc,user->ni,&user->ci);CHKERRQ(ierr); /* a 2x1 vec for inequality constraints */ 161aad13602SShrirang Abhyankar 162aad13602SShrirang Abhyankar /* nexn & nixn matricies for equaly and inequalty constriants */ 163*8e58fa1dSresundermann if (!user->noeqflag){ 164aad13602SShrirang Abhyankar ierr = MatCreate(PETSC_COMM_WORLD,&user->Ae);CHKERRQ(ierr); 165*8e58fa1dSresundermann ierr = MatSetSizes(user->Ae,neloc,nloc,user->ne,user->n);CHKERRQ(ierr); 166*8e58fa1dSresundermann ierr = MatSetUp(user->Ae);CHKERRQ(ierr); 167*8e58fa1dSresundermann ierr = MatSetFromOptions(user->Ae);CHKERRQ(ierr); 168*8e58fa1dSresundermann } 169*8e58fa1dSresundermann 170aad13602SShrirang Abhyankar ierr = MatCreate(PETSC_COMM_WORLD,&user->Ai);CHKERRQ(ierr); 171aad13602SShrirang Abhyankar ierr = MatCreate(PETSC_COMM_WORLD,&user->H);CHKERRQ(ierr); 172aad13602SShrirang Abhyankar 173aad13602SShrirang Abhyankar ierr = MatSetSizes(user->Ai,niloc,nloc,user->ni,user->n);CHKERRQ(ierr); 174aad13602SShrirang Abhyankar ierr = MatSetSizes(user->H,nloc,nloc,user->n,user->n);CHKERRQ(ierr); 175aad13602SShrirang Abhyankar 176aad13602SShrirang Abhyankar ierr = MatSetUp(user->Ai);CHKERRQ(ierr); 177aad13602SShrirang Abhyankar ierr = MatSetUp(user->H);CHKERRQ(ierr); 178aad13602SShrirang Abhyankar 179aad13602SShrirang Abhyankar ierr = MatSetFromOptions(user->Ai);CHKERRQ(ierr); 180aad13602SShrirang Abhyankar ierr = MatSetFromOptions(user->H);CHKERRQ(ierr); 181aad13602SShrirang Abhyankar PetscFunctionReturn(0); 182aad13602SShrirang Abhyankar } 183aad13602SShrirang Abhyankar 184aad13602SShrirang Abhyankar PetscErrorCode DestroyProblem(AppCtx *user) 185aad13602SShrirang Abhyankar { 186aad13602SShrirang Abhyankar PetscErrorCode ierr; 187aad13602SShrirang Abhyankar 188aad13602SShrirang Abhyankar PetscFunctionBegin; 189*8e58fa1dSresundermann if (!user->noeqflag){ 190aad13602SShrirang Abhyankar ierr = MatDestroy(&user->Ae);CHKERRQ(ierr); 191*8e58fa1dSresundermann } 192aad13602SShrirang Abhyankar ierr = MatDestroy(&user->Ai);CHKERRQ(ierr); 193aad13602SShrirang Abhyankar ierr = MatDestroy(&user->H);CHKERRQ(ierr); 194aad13602SShrirang Abhyankar 195aad13602SShrirang Abhyankar ierr = VecDestroy(&user->x);CHKERRQ(ierr); 196*8e58fa1dSresundermann if (!user->noeqflag){ 197aad13602SShrirang Abhyankar ierr = VecDestroy(&user->ce);CHKERRQ(ierr); 198*8e58fa1dSresundermann } 199aad13602SShrirang Abhyankar ierr = VecDestroy(&user->ci);CHKERRQ(ierr); 200aad13602SShrirang Abhyankar ierr = VecDestroy(&user->xl);CHKERRQ(ierr); 201aad13602SShrirang Abhyankar ierr = VecDestroy(&user->xu);CHKERRQ(ierr); 202aad13602SShrirang Abhyankar ierr = VecDestroy(&user->Xseq);CHKERRQ(ierr); 203aad13602SShrirang Abhyankar ierr = VecScatterDestroy(&user->scat);CHKERRQ(ierr); 204aad13602SShrirang Abhyankar PetscFunctionReturn(0); 205aad13602SShrirang Abhyankar } 206aad13602SShrirang Abhyankar 207aad13602SShrirang Abhyankar /* 208aad13602SShrirang Abhyankar f(X) = (x0 - 2)^2 + (x1 - 2)^2 - 2*(x0 + x1) 209aad13602SShrirang Abhyankar G(X) = fx = [2.0*(x[0]-2.0) - 2.0; 2.0*(x[1]-2.0) - 2.0] 210aad13602SShrirang Abhyankar */ 211aad13602SShrirang Abhyankar PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ctx) 212aad13602SShrirang Abhyankar { 213aad13602SShrirang Abhyankar PetscScalar g; 214aad13602SShrirang Abhyankar const PetscScalar *x; 215aad13602SShrirang Abhyankar MPI_Comm comm; 216aad13602SShrirang Abhyankar PetscMPIInt rank; 217aad13602SShrirang Abhyankar PetscErrorCode ierr; 218aad13602SShrirang Abhyankar PetscReal fin; 219aad13602SShrirang Abhyankar AppCtx *user=(AppCtx*)ctx; 220aad13602SShrirang Abhyankar Vec Xseq=user->Xseq; 221aad13602SShrirang Abhyankar VecScatter scat=user->scat; 222aad13602SShrirang Abhyankar 223aad13602SShrirang Abhyankar PetscFunctionBegin; 224aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 225aad13602SShrirang Abhyankar ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 226aad13602SShrirang Abhyankar 227aad13602SShrirang Abhyankar ierr = VecScatterBegin(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 228aad13602SShrirang Abhyankar ierr = VecScatterEnd(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 229aad13602SShrirang Abhyankar 230aad13602SShrirang Abhyankar fin = 0.0; 231aad13602SShrirang Abhyankar if (!rank) { 232aad13602SShrirang Abhyankar ierr = VecGetArrayRead(Xseq,&x);CHKERRQ(ierr); 233aad13602SShrirang Abhyankar fin = (x[0]-2.0)*(x[0]-2.0) + (x[1]-2.0)*(x[1]-2.0) - 2.0*(x[0]+x[1]); 234aad13602SShrirang Abhyankar g = 2.0*(x[0]-2.0) - 2.0; 235aad13602SShrirang Abhyankar ierr = VecSetValue(G,0,g,INSERT_VALUES);CHKERRQ(ierr); 236aad13602SShrirang Abhyankar g = 2.0*(x[1]-2.0) - 2.0; 237aad13602SShrirang Abhyankar ierr = VecSetValue(G,1,g,INSERT_VALUES);CHKERRQ(ierr); 238aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(Xseq,&x);CHKERRQ(ierr); 239aad13602SShrirang Abhyankar } 240aad13602SShrirang Abhyankar ierr = MPI_Allreduce(&fin,f,1,MPIU_REAL,MPIU_SUM,comm);CHKERRQ(ierr); 241aad13602SShrirang Abhyankar ierr = VecAssemblyBegin(G);CHKERRQ(ierr); 242aad13602SShrirang Abhyankar ierr = VecAssemblyEnd(G);CHKERRQ(ierr); 243aad13602SShrirang Abhyankar PetscFunctionReturn(0); 244aad13602SShrirang Abhyankar } 245aad13602SShrirang Abhyankar 246aad13602SShrirang Abhyankar PetscErrorCode FormHessian(Tao tao, Vec x,Mat H, Mat Hpre, void *ctx) 247aad13602SShrirang Abhyankar { 248aad13602SShrirang Abhyankar AppCtx *user=(AppCtx*)ctx; 249aad13602SShrirang Abhyankar Vec DE,DI; 250aad13602SShrirang Abhyankar const PetscScalar *de, *di; 251aad13602SShrirang Abhyankar PetscInt zero=0,one=1; 252aad13602SShrirang Abhyankar PetscScalar two=2.0; 253aad13602SShrirang Abhyankar PetscScalar val=0.0; 254aad13602SShrirang Abhyankar Vec Deseq,Diseq=user->Xseq; 255aad13602SShrirang Abhyankar VecScatter Descat,Discat=user->scat; 256aad13602SShrirang Abhyankar PetscMPIInt rank; 257aad13602SShrirang Abhyankar MPI_Comm comm; 258aad13602SShrirang Abhyankar PetscErrorCode ierr; 259aad13602SShrirang Abhyankar 260aad13602SShrirang Abhyankar PetscFunctionBegin; 261aad13602SShrirang Abhyankar ierr = TaoGetDualVariables(tao,&DE,&DI);CHKERRQ(ierr); 262aad13602SShrirang Abhyankar 263aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 264aad13602SShrirang Abhyankar ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 265aad13602SShrirang Abhyankar 266*8e58fa1dSresundermann if (!user->noeqflag){ 267aad13602SShrirang Abhyankar ierr = VecScatterCreateToZero(DE,&Descat,&Deseq);CHKERRQ(ierr); 268aad13602SShrirang Abhyankar ierr = VecScatterBegin(Descat,DE,Deseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 269aad13602SShrirang Abhyankar ierr = VecScatterEnd(Descat,DE,Deseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 270*8e58fa1dSresundermann } 271*8e58fa1dSresundermann 272*8e58fa1dSresundermann ierr = VecScatterBegin(Discat,DI,Diseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 273aad13602SShrirang Abhyankar ierr = VecScatterEnd(Discat,DI,Diseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 274aad13602SShrirang Abhyankar 275*8e58fa1dSresundermann 276aad13602SShrirang Abhyankar if (!rank){ 277*8e58fa1dSresundermann if (!user->noeqflag){ 278aad13602SShrirang Abhyankar ierr = VecGetArrayRead(Deseq,&de);CHKERRQ(ierr); /* places equality constraint dual into array */ 279*8e58fa1dSresundermann } 280*8e58fa1dSresundermann 281aad13602SShrirang Abhyankar ierr = VecGetArrayRead(Diseq,&di);CHKERRQ(ierr); /* places inequality constraint dual into array */ 282*8e58fa1dSresundermann 283*8e58fa1dSresundermann if (!user->noeqflag){ 284aad13602SShrirang Abhyankar val = 2.0 * (1 + de[0] + di[0] - di[1]); 285aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(Deseq,&de);CHKERRQ(ierr); 286*8e58fa1dSresundermann }else{ 287*8e58fa1dSresundermann val = 2.0 * (1 + di[0] - di[1]); 288*8e58fa1dSresundermann } 289aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(Diseq,&di);CHKERRQ(ierr); 290aad13602SShrirang Abhyankar ierr = MatSetValues(H,1,&zero,1,&zero,&val,INSERT_VALUES);CHKERRQ(ierr); 291aad13602SShrirang Abhyankar ierr = MatSetValues(H,1,&one,1,&one,&two,INSERT_VALUES);CHKERRQ(ierr); 292aad13602SShrirang Abhyankar } 293aad13602SShrirang Abhyankar 294aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 295aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 296*8e58fa1dSresundermann if (!user->noeqflag){ 297aad13602SShrirang Abhyankar ierr = VecScatterDestroy(&Descat);CHKERRQ(ierr); 298aad13602SShrirang Abhyankar ierr = VecDestroy(&Deseq);CHKERRQ(ierr); 299*8e58fa1dSresundermann } 300aad13602SShrirang Abhyankar PetscFunctionReturn(0); 301aad13602SShrirang Abhyankar } 302aad13602SShrirang Abhyankar 303aad13602SShrirang Abhyankar PetscErrorCode FormInequalityConstraints(Tao tao,Vec X,Vec CI,void *ctx) 304aad13602SShrirang Abhyankar { 305aad13602SShrirang Abhyankar const PetscScalar *x; 306aad13602SShrirang Abhyankar PetscScalar ci; 307aad13602SShrirang Abhyankar PetscErrorCode ierr; 308aad13602SShrirang Abhyankar MPI_Comm comm; 309aad13602SShrirang Abhyankar PetscMPIInt rank; 310aad13602SShrirang Abhyankar AppCtx *user=(AppCtx*)ctx; 311aad13602SShrirang Abhyankar Vec Xseq=user->Xseq; 312aad13602SShrirang Abhyankar VecScatter scat=user->scat; 313aad13602SShrirang Abhyankar 314aad13602SShrirang Abhyankar PetscFunctionBegin; 315aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 316aad13602SShrirang Abhyankar ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 317aad13602SShrirang Abhyankar 318aad13602SShrirang Abhyankar ierr = VecScatterBegin(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 319aad13602SShrirang Abhyankar ierr = VecScatterEnd(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 320aad13602SShrirang Abhyankar 321aad13602SShrirang Abhyankar if (!rank) { 322aad13602SShrirang Abhyankar ierr = VecGetArrayRead(Xseq,&x);CHKERRQ(ierr); 323aad13602SShrirang Abhyankar ci = x[0]*x[0] - x[1]; 324aad13602SShrirang Abhyankar ierr = VecSetValue(CI,0,ci,INSERT_VALUES);CHKERRQ(ierr); 325aad13602SShrirang Abhyankar ci = -x[0]*x[0] + x[1] + 1.0; 326aad13602SShrirang Abhyankar ierr = VecSetValue(CI,1,ci,INSERT_VALUES);CHKERRQ(ierr); 327aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(Xseq,&x);CHKERRQ(ierr); 328aad13602SShrirang Abhyankar } 329aad13602SShrirang Abhyankar ierr = VecAssemblyBegin(CI);CHKERRQ(ierr); 330aad13602SShrirang Abhyankar ierr = VecAssemblyEnd(CI);CHKERRQ(ierr); 331aad13602SShrirang Abhyankar PetscFunctionReturn(0); 332aad13602SShrirang Abhyankar } 333aad13602SShrirang Abhyankar 334aad13602SShrirang Abhyankar PetscErrorCode FormEqualityConstraints(Tao tao,Vec X,Vec CE,void *ctx) 335aad13602SShrirang Abhyankar { 336aad13602SShrirang Abhyankar const PetscScalar *x; 337aad13602SShrirang Abhyankar PetscScalar ce; 338aad13602SShrirang Abhyankar PetscErrorCode ierr; 339aad13602SShrirang Abhyankar MPI_Comm comm; 340aad13602SShrirang Abhyankar PetscMPIInt rank; 341aad13602SShrirang Abhyankar AppCtx *user=(AppCtx*)ctx; 342aad13602SShrirang Abhyankar Vec Xseq=user->Xseq; 343aad13602SShrirang Abhyankar VecScatter scat=user->scat; 344aad13602SShrirang Abhyankar 345aad13602SShrirang Abhyankar PetscFunctionBegin; 346aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 347aad13602SShrirang Abhyankar ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 348aad13602SShrirang Abhyankar 349aad13602SShrirang Abhyankar ierr = VecScatterBegin(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 350aad13602SShrirang Abhyankar ierr = VecScatterEnd(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 351aad13602SShrirang Abhyankar 352aad13602SShrirang Abhyankar if (!rank) { 353aad13602SShrirang Abhyankar ierr = VecGetArrayRead(Xseq,&x);CHKERRQ(ierr); 354aad13602SShrirang Abhyankar ce = x[0]*x[0] + x[1] - 2.0; 355aad13602SShrirang Abhyankar ierr = VecSetValue(CE,0,ce,INSERT_VALUES);CHKERRQ(ierr); 356aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(Xseq,&x);CHKERRQ(ierr); 357aad13602SShrirang Abhyankar } 358aad13602SShrirang Abhyankar ierr = VecAssemblyBegin(CE);CHKERRQ(ierr); 359aad13602SShrirang Abhyankar ierr = VecAssemblyEnd(CE);CHKERRQ(ierr); 360aad13602SShrirang Abhyankar PetscFunctionReturn(0); 361aad13602SShrirang Abhyankar } 362aad13602SShrirang Abhyankar 363aad13602SShrirang Abhyankar PetscErrorCode FormInequalityJacobian(Tao tao, Vec X, Mat JI, Mat JIpre, void *ctx) 364aad13602SShrirang Abhyankar { 365aad13602SShrirang Abhyankar AppCtx *user=(AppCtx*)ctx; 366aad13602SShrirang Abhyankar PetscInt cols[2],min,max,i; 367aad13602SShrirang Abhyankar PetscScalar vals[2]; 368aad13602SShrirang Abhyankar const PetscScalar *x; 369aad13602SShrirang Abhyankar PetscErrorCode ierr; 370aad13602SShrirang Abhyankar Vec Xseq=user->Xseq; 371aad13602SShrirang Abhyankar VecScatter scat=user->scat; 372aad13602SShrirang Abhyankar MPI_Comm comm; 373aad13602SShrirang Abhyankar PetscMPIInt rank; 374aad13602SShrirang Abhyankar 375aad13602SShrirang Abhyankar PetscFunctionBegin; 376aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 377aad13602SShrirang Abhyankar ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 378aad13602SShrirang Abhyankar ierr = VecScatterBegin(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 379aad13602SShrirang Abhyankar ierr = VecScatterEnd(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 380aad13602SShrirang Abhyankar 381aad13602SShrirang Abhyankar ierr = VecGetArrayRead(Xseq,&x);CHKERRQ(ierr); 382aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(JI,&min,&max);CHKERRQ(ierr); 383aad13602SShrirang Abhyankar 384aad13602SShrirang Abhyankar cols[0] = 0; cols[1] = 1; 385aad13602SShrirang Abhyankar for (i=min;i<max;i++) { 386aad13602SShrirang Abhyankar if (i==0){ 387aad13602SShrirang Abhyankar vals[0] = +2*x[0]; vals[1] = -1.0; 388aad13602SShrirang Abhyankar ierr = MatSetValues(JI,1,&i,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 389aad13602SShrirang Abhyankar } 390aad13602SShrirang Abhyankar if (i==1) { 391aad13602SShrirang Abhyankar vals[0] = -2*x[0]; vals[1] = +1.0; 392aad13602SShrirang Abhyankar ierr = MatSetValues(JI,1,&i,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 393aad13602SShrirang Abhyankar } 394aad13602SShrirang Abhyankar } 395aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(Xseq,&x);CHKERRQ(ierr); 396aad13602SShrirang Abhyankar 397aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(JI,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 398aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(JI,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 399aad13602SShrirang Abhyankar PetscFunctionReturn(0); 400aad13602SShrirang Abhyankar } 401aad13602SShrirang Abhyankar 402aad13602SShrirang Abhyankar PetscErrorCode FormEqualityJacobian(Tao tao,Vec X,Mat JE,Mat JEpre,void *ctx) 403aad13602SShrirang Abhyankar { 404aad13602SShrirang Abhyankar PetscInt rows[2]; 405aad13602SShrirang Abhyankar PetscScalar vals[2]; 406aad13602SShrirang Abhyankar const PetscScalar *x; 407aad13602SShrirang Abhyankar PetscMPIInt rank; 408aad13602SShrirang Abhyankar MPI_Comm comm; 409aad13602SShrirang Abhyankar PetscErrorCode ierr; 410aad13602SShrirang Abhyankar 411aad13602SShrirang Abhyankar PetscFunctionBegin; 412aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 413aad13602SShrirang Abhyankar ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 414aad13602SShrirang Abhyankar 415aad13602SShrirang Abhyankar if (!rank) { 416aad13602SShrirang Abhyankar ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 417aad13602SShrirang Abhyankar rows[0] = 0; rows[1] = 1; 418aad13602SShrirang Abhyankar vals[0] = 2*x[0]; vals[1] = 1.0; 419aad13602SShrirang Abhyankar ierr = MatSetValues(JE,1,rows,2,rows,vals,INSERT_VALUES);CHKERRQ(ierr); 420aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 421aad13602SShrirang Abhyankar } 422aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(JE,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 423aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(JE,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 424aad13602SShrirang Abhyankar PetscFunctionReturn(0); 425aad13602SShrirang Abhyankar } 426aad13602SShrirang Abhyankar 427aad13602SShrirang Abhyankar 428aad13602SShrirang Abhyankar /*TEST 429aad13602SShrirang Abhyankar 430aad13602SShrirang Abhyankar build: 431aad13602SShrirang Abhyankar requires: !complex !define(PETSC_USE_CXX) 432aad13602SShrirang Abhyankar 433aad13602SShrirang Abhyankar test: 434aad13602SShrirang Abhyankar requires: superlu_dist 435aad13602SShrirang Abhyankar args: -tao_converged_reason 436aad13602SShrirang Abhyankar 437aad13602SShrirang Abhyankar test: 438aad13602SShrirang Abhyankar suffix: 2 439aad13602SShrirang Abhyankar requires: superlu_dist 440aad13602SShrirang Abhyankar nsize: 2 441aad13602SShrirang Abhyankar args: -tao_converged_reason 442aad13602SShrirang Abhyankar 443*8e58fa1dSresundermann test: 444*8e58fa1dSresundermann suffix: 3 445*8e58fa1dSresundermann requires: superlu_dist 446*8e58fa1dSresundermann args: -tao_converged_reason -no_eq 447*8e58fa1dSresundermann 448*8e58fa1dSresundermann test: 449*8e58fa1dSresundermann suffix: 4 450*8e58fa1dSresundermann requires: superlu_dist 451*8e58fa1dSresundermann nsize: 2 452*8e58fa1dSresundermann args: -tao_converged_reason -no_eq 453*8e58fa1dSresundermann 454aad13602SShrirang Abhyankar TEST*/ 455