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