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