1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Finds optimal parameter P_m for the generator system while maintaining generator stability.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*F 5c4762a1bSJed Brown 6c4762a1bSJed Brown \begin{eqnarray} 7c4762a1bSJed Brown \frac{d \theta}{dt} = \omega_b (\omega - \omega_s) 8c4762a1bSJed Brown \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\ 9c4762a1bSJed Brown \end{eqnarray} 10c4762a1bSJed Brown 11c4762a1bSJed Brown F*/ 12c4762a1bSJed Brown 13c4762a1bSJed Brown /* 14c4762a1bSJed Brown Solve the same optimization problem as in ex3opt.c. 15c4762a1bSJed Brown Use finite difference to approximate the gradients. 16c4762a1bSJed Brown */ 17c4762a1bSJed Brown #include <petsctao.h> 18c4762a1bSJed Brown #include <petscts.h> 19c4762a1bSJed Brown #include "ex3.h" 20c4762a1bSJed Brown 21c4762a1bSJed Brown PetscErrorCode FormFunction(Tao,Vec,PetscReal*,void*); 22c4762a1bSJed Brown 23c4762a1bSJed Brown PetscErrorCode monitor(Tao tao,AppCtx *ctx) 24c4762a1bSJed Brown { 25c4762a1bSJed Brown FILE *fp; 26c4762a1bSJed Brown PetscInt iterate; 27c4762a1bSJed Brown PetscReal f,gnorm,cnorm,xdiff; 28c4762a1bSJed Brown Vec X,G; 29c4762a1bSJed Brown const PetscScalar *x,*g; 30c4762a1bSJed Brown TaoConvergedReason reason; 31c4762a1bSJed Brown PetscErrorCode ierr; 32c4762a1bSJed Brown 33c4762a1bSJed Brown PetscFunctionBeginUser; 34c4762a1bSJed Brown ierr = TaoGetSolutionStatus(tao,&iterate,&f,&gnorm,&cnorm,&xdiff,&reason);CHKERRQ(ierr); 35a82e8c82SStefano Zampini ierr = TaoGetSolution(tao,&X);CHKERRQ(ierr); 36a82e8c82SStefano Zampini ierr = TaoGetGradient(tao,&G,NULL,NULL);CHKERRQ(ierr); 37c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = VecGetArrayRead(G,&g);CHKERRQ(ierr); 39c4762a1bSJed Brown fp = fopen("ex3opt_fd_conv.out","a"); 40c4762a1bSJed Brown ierr = PetscFPrintf(PETSC_COMM_WORLD,fp,"%d %g %.12lf %.12lf\n",iterate,gnorm,x[0],g[0]);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = VecRestoreArrayRead(G,&g);CHKERRQ(ierr); 43c4762a1bSJed Brown fclose(fp); 44c4762a1bSJed Brown PetscFunctionReturn(0); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown int main(int argc,char **argv) 48c4762a1bSJed Brown { 49c4762a1bSJed Brown Vec p; 50c4762a1bSJed Brown PetscScalar *x_ptr; 51c4762a1bSJed Brown PetscErrorCode ierr; 52c4762a1bSJed Brown PetscMPIInt size; 53c4762a1bSJed Brown AppCtx ctx; 54c4762a1bSJed Brown Vec lowerb,upperb; 55c4762a1bSJed Brown Tao tao; 56c4762a1bSJed Brown KSP ksp; 57c4762a1bSJed Brown PC pc; 58c4762a1bSJed Brown PetscBool printtofile; 59c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 60c4762a1bSJed Brown Initialize program 61c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 62c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 63c4762a1bSJed Brown PetscFunctionBeginUser; 64ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 65*3c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 68c4762a1bSJed Brown Set runtime options 69c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 70c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr); 71c4762a1bSJed Brown { 72c4762a1bSJed Brown ctx.beta = 2; 73c4762a1bSJed Brown ctx.c = 10000.0; 74c4762a1bSJed Brown ctx.u_s = 1.0; 75c4762a1bSJed Brown ctx.omega_s = 1.0; 76c4762a1bSJed Brown ctx.omega_b = 120.0*PETSC_PI; 77c4762a1bSJed Brown ctx.H = 5.0; 78c4762a1bSJed Brown ierr = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr); 79c4762a1bSJed Brown ctx.D = 5.0; 80c4762a1bSJed Brown ierr = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr); 81c4762a1bSJed Brown ctx.E = 1.1378; 82c4762a1bSJed Brown ctx.V = 1.0; 83c4762a1bSJed Brown ctx.X = 0.545; 84c4762a1bSJed Brown ctx.Pmax = ctx.E*ctx.V/ctx.X; 85c4762a1bSJed Brown ctx.Pmax_ini = ctx.Pmax; 86c4762a1bSJed Brown ierr = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr); 87c4762a1bSJed Brown ctx.Pm = 1.06; 88c4762a1bSJed Brown ierr = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr); 89c4762a1bSJed Brown ctx.tf = 0.1; 90c4762a1bSJed Brown ctx.tcl = 0.2; 91c4762a1bSJed Brown ierr = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr); 93c4762a1bSJed Brown printtofile = PETSC_FALSE; 94c4762a1bSJed Brown ierr = PetscOptionsBool("-printtofile","Print convergence results to file","",printtofile,&printtofile,NULL);CHKERRQ(ierr); 95c4762a1bSJed Brown } 96c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 99c4762a1bSJed Brown ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr); 101c4762a1bSJed Brown if (printtofile) { 102c4762a1bSJed Brown ierr = TaoSetMonitor(tao,(PetscErrorCode (*)(Tao, void*))monitor,(void *)&ctx,PETSC_NULL);CHKERRQ(ierr); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown ierr = TaoSetMaximumIterations(tao,30);CHKERRQ(ierr); 105c4762a1bSJed Brown /* 106c4762a1bSJed Brown Optimization starts 107c4762a1bSJed Brown */ 108c4762a1bSJed Brown /* Set initial solution guess */ 109c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_WORLD,1,&p);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = VecGetArray(p,&x_ptr);CHKERRQ(ierr); 111c4762a1bSJed Brown x_ptr[0] = ctx.Pm; 112c4762a1bSJed Brown ierr = VecRestoreArray(p,&x_ptr);CHKERRQ(ierr); 113c4762a1bSJed Brown 114a82e8c82SStefano Zampini ierr = TaoSetSolution(tao,p);CHKERRQ(ierr); 115c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 116a82e8c82SStefano Zampini ierr = TaoSetObjective(tao,FormFunction,(void *)&ctx);CHKERRQ(ierr); 117a82e8c82SStefano Zampini ierr = TaoSetGradient(tao,NULL,TaoDefaultComputeGradient,(void *)&ctx);CHKERRQ(ierr); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* Set bounds for the optimization */ 120c4762a1bSJed Brown ierr = VecDuplicate(p,&lowerb);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = VecDuplicate(p,&upperb);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = VecGetArray(lowerb,&x_ptr);CHKERRQ(ierr); 123c4762a1bSJed Brown x_ptr[0] = 0.; 124c4762a1bSJed Brown ierr = VecRestoreArray(lowerb,&x_ptr);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = VecGetArray(upperb,&x_ptr);CHKERRQ(ierr); 126c4762a1bSJed Brown x_ptr[0] = 1.1; 127c4762a1bSJed Brown ierr = VecRestoreArray(upperb,&x_ptr);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = TaoSetVariableBounds(tao,lowerb,upperb);CHKERRQ(ierr); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* Check for any TAO command line options */ 131c4762a1bSJed Brown ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr); 133c4762a1bSJed Brown if (ksp) { 134c4762a1bSJed Brown ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 139c4762a1bSJed Brown ierr = TaoSolve(tao);CHKERRQ(ierr); 140c4762a1bSJed Brown 141c4762a1bSJed Brown ierr = VecView(p,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* Free TAO data structures */ 144c4762a1bSJed Brown ierr = TaoDestroy(&tao);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = VecDestroy(&p);CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = VecDestroy(&lowerb);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = VecDestroy(&upperb);CHKERRQ(ierr); 148c4762a1bSJed Brown ierr = PetscFinalize(); 149c4762a1bSJed Brown return ierr; 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 153c4762a1bSJed Brown /* 154c4762a1bSJed Brown FormFunction - Evaluates the function and corresponding gradient. 155c4762a1bSJed Brown 156c4762a1bSJed Brown Input Parameters: 157c4762a1bSJed Brown tao - the Tao context 158c4762a1bSJed Brown X - the input vector 159a82e8c82SStefano Zampini ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient() 160c4762a1bSJed Brown 161c4762a1bSJed Brown Output Parameters: 162c4762a1bSJed Brown f - the newly evaluated function 163c4762a1bSJed Brown */ 164c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0) 165c4762a1bSJed Brown { 166c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)ctx0; 167c4762a1bSJed Brown TS ts,quadts; 168c4762a1bSJed Brown Vec U; /* solution will be stored here */ 169c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 170c4762a1bSJed Brown PetscErrorCode ierr; 171c4762a1bSJed Brown PetscInt n = 2; 172c4762a1bSJed Brown PetscReal ftime; 173c4762a1bSJed Brown PetscInt steps; 174c4762a1bSJed Brown PetscScalar *u; 175c4762a1bSJed Brown const PetscScalar *x_ptr,*qx_ptr; 176c4762a1bSJed Brown Vec q; 177c4762a1bSJed Brown PetscInt direction[2]; 178c4762a1bSJed Brown PetscBool terminate[2]; 179c4762a1bSJed Brown 180c4762a1bSJed Brown ierr = VecGetArrayRead(P,&x_ptr);CHKERRQ(ierr); 181c4762a1bSJed Brown ctx->Pm = x_ptr[0]; 182c4762a1bSJed Brown ierr = VecRestoreArrayRead(P,&x_ptr);CHKERRQ(ierr); 183c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 184c4762a1bSJed Brown Create necessary matrix and vectors 185c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 186c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 188c4762a1bSJed Brown ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr); 189c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 191c4762a1bSJed Brown 192c4762a1bSJed Brown ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr); 193c4762a1bSJed Brown 194c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 195c4762a1bSJed Brown Create timestepping solver context 196c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 197c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = TSSetType(ts,TSCN);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);CHKERRQ(ierr); 201c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,ctx);CHKERRQ(ierr); 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 204c4762a1bSJed Brown Set initial conditions 205c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 206c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 207c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax); 208c4762a1bSJed Brown u[1] = 1.0; 209c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = TSSetSolution(ts,U);CHKERRQ(ierr); 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 213c4762a1bSJed Brown Set solver options 214c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 215c4762a1bSJed Brown ierr = TSSetMaxTime(ts,1.0);CHKERRQ(ierr); 216c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.03125);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = TSCreateQuadratureTS(ts,PETSC_TRUE,&quadts);CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = TSGetSolution(quadts,&q);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = VecSet(q,0.0);CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,ctx);CHKERRQ(ierr); 222c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 223c4762a1bSJed Brown 224c4762a1bSJed Brown direction[0] = direction[1] = 1; 225c4762a1bSJed Brown terminate[0] = terminate[1] = PETSC_FALSE; 226c4762a1bSJed Brown 227c4762a1bSJed Brown ierr = TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)ctx);CHKERRQ(ierr); 228c4762a1bSJed Brown 229c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 230c4762a1bSJed Brown Solve nonlinear system 231c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 232c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 233c4762a1bSJed Brown 234c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 236c4762a1bSJed Brown ierr = VecGetArrayRead(q,&qx_ptr);CHKERRQ(ierr); 237c4762a1bSJed Brown *f = -ctx->Pm + qx_ptr[0]; 238c4762a1bSJed Brown ierr = VecRestoreArrayRead(q,&qx_ptr);CHKERRQ(ierr); 239c4762a1bSJed Brown 240c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 241c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 242c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 243c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 245c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 246c4762a1bSJed Brown return 0; 247c4762a1bSJed Brown } 248c4762a1bSJed Brown 249c4762a1bSJed Brown /*TEST 250c4762a1bSJed Brown 251c4762a1bSJed Brown build: 252c4762a1bSJed Brown requires: !complex !single 253c4762a1bSJed Brown 254c4762a1bSJed Brown test: 255c4762a1bSJed Brown args: -ts_type cn -pc_type lu -tao_monitor -tao_gatol 1e-3 256c4762a1bSJed Brown 257c4762a1bSJed Brown TEST*/ 258