1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Basic equation for generator stability analysis.\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 Ensemble of initial conditions 12c4762a1bSJed Brown ./ex2 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly 13c4762a1bSJed Brown 14c4762a1bSJed Brown Fault at .1 seconds 15c4762a1bSJed Brown ./ex2 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly 16c4762a1bSJed Brown 17c4762a1bSJed Brown Initial conditions same as when fault is ended 18c4762a1bSJed Brown ./ex2 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly 19c4762a1bSJed Brown 20c4762a1bSJed Brown F*/ 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* 23c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this 24c4762a1bSJed Brown file automatically includes: 25c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 26c4762a1bSJed Brown petscmat.h - matrices 27c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 28c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 29c4762a1bSJed Brown petscksp.h - linear solvers 30c4762a1bSJed Brown */ 31c4762a1bSJed Brown 32c4762a1bSJed Brown #include <petsctao.h> 33c4762a1bSJed Brown #include <petscts.h> 34c4762a1bSJed Brown 35c4762a1bSJed Brown typedef struct { 36c4762a1bSJed Brown TS ts; 37c4762a1bSJed Brown PetscScalar H,D,omega_b,omega_s,Pmax,Pm,E,V,X,u_s,c; 38c4762a1bSJed Brown PetscInt beta; 39c4762a1bSJed Brown PetscReal tf,tcl,dt; 40c4762a1bSJed Brown } AppCtx; 41c4762a1bSJed Brown 42c4762a1bSJed Brown PetscErrorCode FormFunction(Tao,Vec,PetscReal*,void*); 43c4762a1bSJed Brown PetscErrorCode FormGradient(Tao,Vec,Vec,void*); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* 46c4762a1bSJed Brown Defines the ODE passed to the ODE solver 47c4762a1bSJed Brown */ 48c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 49c4762a1bSJed Brown { 50c4762a1bSJed Brown PetscScalar *f,Pmax; 51c4762a1bSJed Brown const PetscScalar *u; 52c4762a1bSJed Brown 53c4762a1bSJed Brown PetscFunctionBegin; 54c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 55*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 56*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 57c4762a1bSJed Brown if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */ 58c4762a1bSJed Brown else Pmax = ctx->Pmax; 59c4762a1bSJed Brown 60c4762a1bSJed Brown f[0] = ctx->omega_b*(u[1] - ctx->omega_s); 61c4762a1bSJed Brown f[1] = (-Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s) + ctx->Pm)*ctx->omega_s/(2.0*ctx->H); 62c4762a1bSJed Brown 63*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 64*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 65c4762a1bSJed Brown PetscFunctionReturn(0); 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* 69c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 70c4762a1bSJed Brown */ 71c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx) 72c4762a1bSJed Brown { 73c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 74c4762a1bSJed Brown PetscScalar J[2][2],Pmax; 75c4762a1bSJed Brown const PetscScalar *u; 76c4762a1bSJed Brown 77c4762a1bSJed Brown PetscFunctionBegin; 78*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 79c4762a1bSJed Brown if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */ 80c4762a1bSJed Brown else Pmax = ctx->Pmax; 81c4762a1bSJed Brown 82c4762a1bSJed Brown J[0][0] = 0; J[0][1] = ctx->omega_b; 83c4762a1bSJed Brown J[1][1] = -ctx->D*ctx->omega_s/(2.0*ctx->H); J[1][0] = -Pmax*PetscCosScalar(u[0])*ctx->omega_s/(2.0*ctx->H); 84c4762a1bSJed Brown 85*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 86*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 87c4762a1bSJed Brown 88*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 89*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 90c4762a1bSJed Brown if (A != B) { 91*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 92*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown PetscFunctionReturn(0); 95c4762a1bSJed Brown } 96c4762a1bSJed Brown 97c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0) 98c4762a1bSJed Brown { 99c4762a1bSJed Brown PetscInt row[] = {0,1},col[]={0}; 100c4762a1bSJed Brown PetscScalar J[2][1]; 101c4762a1bSJed Brown AppCtx *ctx=(AppCtx*)ctx0; 102c4762a1bSJed Brown 103c4762a1bSJed Brown PetscFunctionBeginUser; 104c4762a1bSJed Brown J[0][0] = 0; 105c4762a1bSJed Brown J[1][0] = ctx->omega_s/(2.0*ctx->H); 106*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES)); 107*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 108*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 109c4762a1bSJed Brown PetscFunctionReturn(0); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112c4762a1bSJed Brown static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx) 113c4762a1bSJed Brown { 114c4762a1bSJed Brown PetscScalar *r; 115c4762a1bSJed Brown const PetscScalar *u; 116c4762a1bSJed Brown 117c4762a1bSJed Brown PetscFunctionBegin; 118*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 119*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(R,&r)); 1202f613bf5SBarry Smith r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta); 121*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(R,&r)); 122*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 123c4762a1bSJed Brown PetscFunctionReturn(0); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,AppCtx *ctx) 127c4762a1bSJed Brown { 128c4762a1bSJed Brown PetscScalar ru[1]; 129c4762a1bSJed Brown const PetscScalar *u; 130c4762a1bSJed Brown PetscInt row[] = {0},col[] = {0}; 131c4762a1bSJed Brown 132c4762a1bSJed Brown PetscFunctionBegin; 133*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 1342f613bf5SBarry Smith ru[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1); 135*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 136*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(DRDU,1,row,1,col,ru,INSERT_VALUES)); 137*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY)); 138*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDU,MAT_FINAL_ASSEMBLY)); 139c4762a1bSJed Brown PetscFunctionReturn(0); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142c4762a1bSJed Brown static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,AppCtx *ctx) 143c4762a1bSJed Brown { 144c4762a1bSJed Brown PetscFunctionBegin; 145*9566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(DRDP)); 146*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY)); 147*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDP,MAT_FINAL_ASSEMBLY)); 148c4762a1bSJed Brown PetscFunctionReturn(0); 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx) 152c4762a1bSJed Brown { 153c4762a1bSJed Brown PetscScalar *y,sensip; 154c4762a1bSJed Brown const PetscScalar *x; 155c4762a1bSJed Brown 156c4762a1bSJed Brown PetscFunctionBegin; 157*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(lambda,&x)); 158*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu,&y)); 159c4762a1bSJed Brown sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0]; 160c4762a1bSJed Brown y[0] = sensip; 161*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu,&y)); 162*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(lambda,&x)); 163c4762a1bSJed Brown PetscFunctionReturn(0); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown int main(int argc,char **argv) 167c4762a1bSJed Brown { 168c4762a1bSJed Brown Vec p; 169c4762a1bSJed Brown PetscScalar *x_ptr; 170c4762a1bSJed Brown PetscErrorCode ierr; 171c4762a1bSJed Brown PetscMPIInt size; 172c4762a1bSJed Brown AppCtx ctx; 173c4762a1bSJed Brown Vec lowerb,upperb; 174c4762a1bSJed Brown Tao tao; 175c4762a1bSJed Brown KSP ksp; 176c4762a1bSJed Brown PC pc; 177c4762a1bSJed Brown Vec U,lambda[1],mu[1]; 178c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 179c4762a1bSJed Brown Mat Jacp; /* Jacobian matrix */ 180c4762a1bSJed Brown Mat DRDU,DRDP; 181c4762a1bSJed Brown PetscInt n = 2; 182c4762a1bSJed Brown TS quadts; 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 185c4762a1bSJed Brown Initialize program 186c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 187*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 188c4762a1bSJed Brown PetscFunctionBeginUser; 189*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1903c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 191c4762a1bSJed Brown 192c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 193c4762a1bSJed Brown Set runtime options 194c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 195*9566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");PetscCall(ierr); 196c4762a1bSJed Brown { 197c4762a1bSJed Brown ctx.beta = 2; 198c4762a1bSJed Brown ctx.c = PetscRealConstant(10000.0); 199c4762a1bSJed Brown ctx.u_s = PetscRealConstant(1.0); 200c4762a1bSJed Brown ctx.omega_s = PetscRealConstant(1.0); 201c4762a1bSJed Brown ctx.omega_b = PetscRealConstant(120.0)*PETSC_PI; 202c4762a1bSJed Brown ctx.H = PetscRealConstant(5.0); 203*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL)); 204c4762a1bSJed Brown ctx.D = PetscRealConstant(5.0); 205*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL)); 206c4762a1bSJed Brown ctx.E = PetscRealConstant(1.1378); 207c4762a1bSJed Brown ctx.V = PetscRealConstant(1.0); 208c4762a1bSJed Brown ctx.X = PetscRealConstant(0.545); 209c4762a1bSJed Brown ctx.Pmax = ctx.E*ctx.V/ctx.X; 210*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL)); 211c4762a1bSJed Brown ctx.Pm = PetscRealConstant(1.0194); 212*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL)); 213c4762a1bSJed Brown ctx.tf = PetscRealConstant(0.1); 214c4762a1bSJed Brown ctx.tcl = PetscRealConstant(0.2); 215*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL)); 216*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL)); 217c4762a1bSJed Brown 218c4762a1bSJed Brown } 219*9566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 220c4762a1bSJed Brown 221c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 222c4762a1bSJed Brown Create necessary matrix and vectors 223c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 224*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 225*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 226*9566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATDENSE)); 227*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 228*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 229c4762a1bSJed Brown 230*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&U,NULL)); 231c4762a1bSJed Brown 232*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&Jacp)); 233*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1)); 234*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jacp)); 235*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jacp)); 236*9566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&DRDP)); 237*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(DRDP)); 238*9566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,2,NULL,&DRDU)); 239*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(DRDU)); 240c4762a1bSJed Brown 241c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 242c4762a1bSJed Brown Create timestepping solver context 243c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 244*9566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ctx.ts)); 245*9566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ctx.ts,TS_NONLINEAR)); 246*9566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(ctx.ts,TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 247*9566063dSJacob Faibussowitsch PetscCall(TSSetType(ctx.ts,TSRK)); 248*9566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ctx.ts,NULL,(TSRHSFunction)RHSFunction,&ctx)); 249*9566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ctx.ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx)); 250*9566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ctx.ts,TS_EXACTFINALTIME_MATCHSTEP)); 251c4762a1bSJed Brown 252*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&lambda[0],NULL)); 253*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Jacp,&mu[0],NULL)); 254*9566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ctx.ts,1,lambda,mu)); 255*9566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ctx.ts,Jacp,RHSJacobianP,&ctx)); 256c4762a1bSJed Brown 257c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 258c4762a1bSJed Brown Set solver options 259c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 260*9566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ctx.ts,PetscRealConstant(1.0))); 261*9566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ctx.ts,PetscRealConstant(0.01))); 262*9566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ctx.ts)); 263c4762a1bSJed Brown 264*9566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ctx.ts,&ctx.dt)); /* save the stepsize */ 265c4762a1bSJed Brown 266*9566063dSJacob Faibussowitsch PetscCall(TSCreateQuadratureTS(ctx.ts,PETSC_TRUE,&quadts)); 267*9566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx)); 268*9566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(quadts,DRDU,DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx)); 269*9566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(quadts,DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,&ctx)); 270*9566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ctx.ts,U)); 271c4762a1bSJed Brown 272c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 273*9566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao)); 274*9566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAOBLMVM)); 275c4762a1bSJed Brown 276c4762a1bSJed Brown /* 277c4762a1bSJed Brown Optimization starts 278c4762a1bSJed Brown */ 279c4762a1bSJed Brown /* Set initial solution guess */ 280*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD,1,&p)); 281*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(p,&x_ptr)); 282c4762a1bSJed Brown x_ptr[0] = ctx.Pm; 283*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(p,&x_ptr)); 284c4762a1bSJed Brown 285*9566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,p)); 286c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 287*9566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao,FormFunction,(void *)&ctx)); 288*9566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao,NULL,FormGradient,(void *)&ctx)); 289c4762a1bSJed Brown 290c4762a1bSJed Brown /* Set bounds for the optimization */ 291*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(p,&lowerb)); 292*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(p,&upperb)); 293*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(lowerb,&x_ptr)); 294c4762a1bSJed Brown x_ptr[0] = 0.; 295*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lowerb,&x_ptr)); 296*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(upperb,&x_ptr)); 297c4762a1bSJed Brown x_ptr[0] = PetscRealConstant(1.1); 298*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(upperb,&x_ptr)); 299*9566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(tao,lowerb,upperb)); 300c4762a1bSJed Brown 301c4762a1bSJed Brown /* Check for any TAO command line options */ 302*9566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 303*9566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao,&ksp)); 304c4762a1bSJed Brown if (ksp) { 305*9566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 306*9566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCNONE)); 307c4762a1bSJed Brown } 308c4762a1bSJed Brown 309c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 310*9566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 311c4762a1bSJed Brown 312*9566063dSJacob Faibussowitsch PetscCall(VecView(p,PETSC_VIEWER_STDOUT_WORLD)); 313c4762a1bSJed Brown /* Free TAO data structures */ 314*9566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 315*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&p)); 316*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lowerb)); 317*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&upperb)); 318c4762a1bSJed Brown 319*9566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ctx.ts)); 320*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 321*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 322*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jacp)); 323*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&DRDU)); 324*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&DRDP)); 325*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0])); 326*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[0])); 327*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 328b122ec5aSJacob Faibussowitsch return 0; 329c4762a1bSJed Brown } 330c4762a1bSJed Brown 331c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 332c4762a1bSJed Brown /* 333c4762a1bSJed Brown FormFunction - Evaluates the function 334c4762a1bSJed Brown 335c4762a1bSJed Brown Input Parameters: 336c4762a1bSJed Brown tao - the Tao context 337c4762a1bSJed Brown X - the input vector 338a82e8c82SStefano Zampini ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient() 339c4762a1bSJed Brown 340c4762a1bSJed Brown Output Parameters: 341c4762a1bSJed Brown f - the newly evaluated function 342c4762a1bSJed Brown */ 343c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0) 344c4762a1bSJed Brown { 345c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)ctx0; 346c4762a1bSJed Brown TS ts = ctx->ts; 347c4762a1bSJed Brown Vec U; /* solution will be stored here */ 348c4762a1bSJed Brown PetscScalar *u; 349c4762a1bSJed Brown PetscScalar *x_ptr; 350c4762a1bSJed Brown Vec q; 351c4762a1bSJed Brown 352*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P,(const PetscScalar**)&x_ptr)); 353c4762a1bSJed Brown ctx->Pm = x_ptr[0]; 354*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr)); 355c4762a1bSJed Brown 356c4762a1bSJed Brown /* reset time */ 357*9566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts,0.0)); 358c4762a1bSJed Brown /* reset step counter, this is critical for adjoint solver */ 359*9566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts,0)); 360c4762a1bSJed Brown /* reset step size, the step size becomes negative after TSAdjointSolve */ 361*9566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,ctx->dt)); 362c4762a1bSJed Brown /* reinitialize the integral value */ 363*9566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts,&q)); 364*9566063dSJacob Faibussowitsch PetscCall(VecSet(q,0.0)); 365c4762a1bSJed Brown 366c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 367c4762a1bSJed Brown Set initial conditions 368c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 369*9566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts,&U)); 370*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&u)); 371c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax); 372c4762a1bSJed Brown u[1] = PetscRealConstant(1.0); 373*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U,&u)); 374c4762a1bSJed Brown 375c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 376c4762a1bSJed Brown Solve nonlinear system 377c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 378*9566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,U)); 379*9566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts,&q)); 380*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(q,&x_ptr)); 381c4762a1bSJed Brown *f = -ctx->Pm + x_ptr[0]; 382*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(q,&x_ptr)); 383c4762a1bSJed Brown return 0; 384c4762a1bSJed Brown } 385c4762a1bSJed Brown 386c4762a1bSJed Brown PetscErrorCode FormGradient(Tao tao,Vec P,Vec G,void *ctx0) 387c4762a1bSJed Brown { 388c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)ctx0; 389c4762a1bSJed Brown TS ts = ctx->ts; 390c4762a1bSJed Brown Vec U; /* solution will be stored here */ 391c4762a1bSJed Brown PetscReal ftime; 392c4762a1bSJed Brown PetscInt steps; 393c4762a1bSJed Brown PetscScalar *u; 394c4762a1bSJed Brown PetscScalar *x_ptr,*y_ptr; 395c4762a1bSJed Brown Vec *lambda,q,*mu; 396c4762a1bSJed Brown 397*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P,(const PetscScalar**)&x_ptr)); 398c4762a1bSJed Brown ctx->Pm = x_ptr[0]; 399*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr)); 400c4762a1bSJed Brown 401c4762a1bSJed Brown /* reset time */ 402*9566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts,0.0)); 403c4762a1bSJed Brown /* reset step counter, this is critical for adjoint solver */ 404*9566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts,0)); 405c4762a1bSJed Brown /* reset step size, the step size becomes negative after TSAdjointSolve */ 406*9566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,ctx->dt)); 407c4762a1bSJed Brown /* reinitialize the integral value */ 408*9566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts,&q)); 409*9566063dSJacob Faibussowitsch PetscCall(VecSet(q,0.0)); 410c4762a1bSJed Brown 411c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 412c4762a1bSJed Brown Set initial conditions 413c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 414*9566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts,&U)); 415*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&u)); 416c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax); 417c4762a1bSJed Brown u[1] = PetscRealConstant(1.0); 418*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U,&u)); 419c4762a1bSJed Brown 420f32d6360SSatish Balay /* Set up to save trajectory before TSSetFromOptions() so that TSTrajectory options can be captured */ 421*9566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 422*9566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 423c4762a1bSJed Brown 424c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 425c4762a1bSJed Brown Solve nonlinear system 426c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 427*9566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,U)); 428c4762a1bSJed Brown 429*9566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ftime)); 430*9566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 431c4762a1bSJed Brown 432c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 433c4762a1bSJed Brown Adjoint model starts here 434c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 435*9566063dSJacob Faibussowitsch PetscCall(TSGetCostGradients(ts,NULL,&lambda,&mu)); 436c4762a1bSJed Brown /* Set initial conditions for the adjoint integration */ 437*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0],&y_ptr)); 438c4762a1bSJed Brown y_ptr[0] = 0.0; y_ptr[1] = 0.0; 439*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0],&y_ptr)); 440*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0],&x_ptr)); 441c4762a1bSJed Brown x_ptr[0] = PetscRealConstant(-1.0); 442*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0],&x_ptr)); 443c4762a1bSJed Brown 444*9566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 445*9566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts,&q)); 446*9566063dSJacob Faibussowitsch PetscCall(ComputeSensiP(lambda[0],mu[0],ctx)); 447*9566063dSJacob Faibussowitsch PetscCall(VecCopy(mu[0],G)); 448c4762a1bSJed Brown return 0; 449c4762a1bSJed Brown } 450c4762a1bSJed Brown 451c4762a1bSJed Brown /*TEST 452c4762a1bSJed Brown 453c4762a1bSJed Brown build: 454c4762a1bSJed Brown requires: !complex 455c4762a1bSJed Brown 456c4762a1bSJed Brown test: 457c4762a1bSJed Brown args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason 458c4762a1bSJed Brown 459c4762a1bSJed Brown test: 460c4762a1bSJed Brown suffix: 2 461c4762a1bSJed Brown args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason -tao_test_gradient 462c4762a1bSJed Brown 463c4762a1bSJed Brown TEST*/ 464