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 12c4762a1bSJed Brown 13c4762a1bSJed Brown Ensemble of initial conditions 14c4762a1bSJed Brown ./ex9 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly 15c4762a1bSJed Brown 16c4762a1bSJed Brown Fault at .1 seconds 17c4762a1bSJed Brown ./ex9 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly 18c4762a1bSJed Brown 19c4762a1bSJed Brown Initial conditions same as when fault is ended 20c4762a1bSJed Brown ./ex9 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly 21c4762a1bSJed Brown 22c4762a1bSJed Brown 23c4762a1bSJed Brown F*/ 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* 26c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this 27c4762a1bSJed Brown file automatically includes: 28c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 29c4762a1bSJed Brown petscmat.h - matrices 30c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 31c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 32c4762a1bSJed Brown petscksp.h - linear solvers 33c4762a1bSJed Brown */ 34c4762a1bSJed Brown 35c4762a1bSJed Brown #include <petscts.h> 36c4762a1bSJed Brown 37c4762a1bSJed Brown typedef struct { 38c4762a1bSJed Brown PetscScalar H,D,omega_b,omega_s,Pmax,Pm,E,V,X,u_s,c; 39c4762a1bSJed Brown PetscInt beta; 40c4762a1bSJed Brown PetscReal tf,tcl; 41c4762a1bSJed Brown } AppCtx; 42c4762a1bSJed Brown 43c4762a1bSJed Brown PetscErrorCode PostStepFunction(TS ts) 44c4762a1bSJed Brown { 45c4762a1bSJed Brown PetscErrorCode ierr; 46c4762a1bSJed Brown Vec U; 47c4762a1bSJed Brown PetscReal t; 48c4762a1bSJed Brown const PetscScalar *u; 49c4762a1bSJed Brown 50c4762a1bSJed Brown PetscFunctionBegin; 51c4762a1bSJed Brown ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"delta(%3.2f) = %8.7f\n",(double)t,(double)u[0]);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 56c4762a1bSJed Brown PetscFunctionReturn(0); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* 60c4762a1bSJed Brown Defines the ODE passed to the ODE solver 61c4762a1bSJed Brown */ 62c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 63c4762a1bSJed Brown { 64c4762a1bSJed Brown PetscErrorCode ierr; 65c4762a1bSJed Brown PetscScalar *f,Pmax; 66c4762a1bSJed Brown const PetscScalar *u; 67c4762a1bSJed Brown 68c4762a1bSJed Brown PetscFunctionBegin; 69c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 70c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 72c4762a1bSJed 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 */ 73c4762a1bSJed Brown else Pmax = ctx->Pmax; 74c4762a1bSJed Brown 75c4762a1bSJed Brown f[0] = ctx->omega_b*(u[1] - ctx->omega_s); 76c4762a1bSJed Brown f[1] = (-Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s) + ctx->Pm)*ctx->omega_s/(2.0*ctx->H); 77c4762a1bSJed Brown 78c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 79c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 80c4762a1bSJed Brown PetscFunctionReturn(0); 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* 84c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 85c4762a1bSJed Brown */ 86c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx) 87c4762a1bSJed Brown { 88c4762a1bSJed Brown PetscErrorCode ierr; 89c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 90c4762a1bSJed Brown PetscScalar J[2][2],Pmax; 91c4762a1bSJed Brown const PetscScalar *u; 92c4762a1bSJed Brown 93c4762a1bSJed Brown PetscFunctionBegin; 94c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 95c4762a1bSJed 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 */ 96c4762a1bSJed Brown else Pmax = ctx->Pmax; 97c4762a1bSJed Brown 98c4762a1bSJed Brown J[0][0] = 0; J[0][1] = ctx->omega_b; 99c4762a1bSJed 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); 100c4762a1bSJed Brown 101c4762a1bSJed Brown ierr = MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 102c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 103c4762a1bSJed Brown 104c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 106c4762a1bSJed Brown if (A != B) { 107c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 109c4762a1bSJed Brown } 110c4762a1bSJed Brown PetscFunctionReturn(0); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0) 114c4762a1bSJed Brown { 115c4762a1bSJed Brown PetscErrorCode ierr; 116c4762a1bSJed Brown PetscInt row[] = {0,1},col[]={0}; 117c4762a1bSJed Brown PetscScalar J[2][1]; 118c4762a1bSJed Brown AppCtx *ctx=(AppCtx*)ctx0; 119c4762a1bSJed Brown 120c4762a1bSJed Brown PetscFunctionBeginUser; 121c4762a1bSJed Brown J[0][0] = 0; 122c4762a1bSJed Brown J[1][0] = ctx->omega_s/(2.0*ctx->H); 123c4762a1bSJed Brown ierr = MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126c4762a1bSJed Brown PetscFunctionReturn(0); 127c4762a1bSJed Brown } 128c4762a1bSJed Brown 129c4762a1bSJed Brown static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx) 130c4762a1bSJed Brown { 131c4762a1bSJed Brown PetscErrorCode ierr; 132c4762a1bSJed Brown PetscScalar *r; 133c4762a1bSJed Brown const PetscScalar *u; 134c4762a1bSJed Brown 135c4762a1bSJed Brown PetscFunctionBegin; 136c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = VecGetArray(R,&r);CHKERRQ(ierr); 138c4762a1bSJed Brown r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = VecRestoreArray(R,&r);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 141c4762a1bSJed Brown PetscFunctionReturn(0); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144c4762a1bSJed Brown static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,AppCtx *ctx) 145c4762a1bSJed Brown { 146c4762a1bSJed Brown PetscErrorCode ierr; 147c4762a1bSJed Brown PetscScalar ru[1]; 148c4762a1bSJed Brown const PetscScalar *u; 149c4762a1bSJed Brown PetscInt row[] = {0},col[] = {0}; 150c4762a1bSJed Brown 151c4762a1bSJed Brown PetscFunctionBegin; 152c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 153c4762a1bSJed Brown ru[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = MatSetValues(DRDU,1,row,1,col,ru,INSERT_VALUES);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = MatAssemblyEnd(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 158c4762a1bSJed Brown PetscFunctionReturn(0); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,AppCtx *ctx) 162c4762a1bSJed Brown { 163c4762a1bSJed Brown PetscErrorCode ierr; 164c4762a1bSJed Brown 165c4762a1bSJed Brown PetscFunctionBegin; 166c4762a1bSJed Brown ierr = MatZeroEntries(DRDP);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = MatAssemblyEnd(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169c4762a1bSJed Brown PetscFunctionReturn(0); 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 172c4762a1bSJed Brown PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx) 173c4762a1bSJed Brown { 174c4762a1bSJed Brown PetscErrorCode ierr; 175c4762a1bSJed Brown PetscScalar sensip; 176c4762a1bSJed Brown const PetscScalar *x,*y; 177c4762a1bSJed Brown 178c4762a1bSJed Brown PetscFunctionBegin; 179c4762a1bSJed Brown ierr = VecGetArrayRead(lambda,&x);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = VecGetArrayRead(mu,&y);CHKERRQ(ierr); 181c4762a1bSJed Brown sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0]; 182c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %.7f \n",(double)sensip);CHKERRQ(ierr); 183c4762a1bSJed Brown ierr = VecRestoreArrayRead(lambda,&x);CHKERRQ(ierr); 184c4762a1bSJed Brown ierr = VecRestoreArrayRead(mu,&y);CHKERRQ(ierr); 185c4762a1bSJed Brown PetscFunctionReturn(0); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 188c4762a1bSJed Brown int main(int argc,char **argv) 189c4762a1bSJed Brown { 190c4762a1bSJed Brown TS ts,quadts; /* ODE integrator */ 191c4762a1bSJed Brown Vec U; /* solution will be stored here */ 192c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 193c4762a1bSJed Brown Mat Jacp; /* Jacobian matrix */ 194c4762a1bSJed Brown Mat DRDU,DRDP; 195c4762a1bSJed Brown PetscErrorCode ierr; 196c4762a1bSJed Brown PetscMPIInt size; 197c4762a1bSJed Brown PetscInt n = 2; 198c4762a1bSJed Brown AppCtx ctx; 199c4762a1bSJed Brown PetscScalar *u; 200c4762a1bSJed Brown PetscReal du[2] = {0.0,0.0}; 201c4762a1bSJed Brown PetscBool ensemble = PETSC_FALSE,flg1,flg2; 202c4762a1bSJed Brown PetscReal ftime; 203c4762a1bSJed Brown PetscInt steps; 204c4762a1bSJed Brown PetscScalar *x_ptr,*y_ptr; 205c4762a1bSJed Brown Vec lambda[1],q,mu[1]; 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 208c4762a1bSJed Brown Initialize program 209c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 210c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 211*ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 212c4762a1bSJed Brown if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs"); 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 215c4762a1bSJed Brown Create necessary matrix and vectors 216c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 217c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 222c4762a1bSJed Brown 223c4762a1bSJed Brown ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr); 224c4762a1bSJed Brown 225c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&Jacp);CHKERRQ(ierr); 226c4762a1bSJed Brown ierr = MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr); 227c4762a1bSJed Brown ierr = MatSetFromOptions(Jacp);CHKERRQ(ierr); 228c4762a1bSJed Brown ierr = MatSetUp(Jacp);CHKERRQ(ierr); 229c4762a1bSJed Brown 230c4762a1bSJed Brown ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&DRDP);CHKERRQ(ierr); 231c4762a1bSJed Brown ierr = MatSetUp(DRDP);CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,2,NULL,&DRDU);CHKERRQ(ierr); 233c4762a1bSJed Brown ierr = MatSetUp(DRDU);CHKERRQ(ierr); 234c4762a1bSJed Brown 235c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 236c4762a1bSJed Brown Set runtime options 237c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 238c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr); 239c4762a1bSJed Brown { 240c4762a1bSJed Brown ctx.beta = 2; 241c4762a1bSJed Brown ctx.c = 10000.0; 242c4762a1bSJed Brown ctx.u_s = 1.0; 243c4762a1bSJed Brown ctx.omega_s = 1.0; 244c4762a1bSJed Brown ctx.omega_b = 120.0*PETSC_PI; 245c4762a1bSJed Brown ctx.H = 5.0; 246c4762a1bSJed Brown ierr = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr); 247c4762a1bSJed Brown ctx.D = 5.0; 248c4762a1bSJed Brown ierr = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr); 249c4762a1bSJed Brown ctx.E = 1.1378; 250c4762a1bSJed Brown ctx.V = 1.0; 251c4762a1bSJed Brown ctx.X = 0.545; 252c4762a1bSJed Brown ctx.Pmax = ctx.E*ctx.V/ctx.X; 253c4762a1bSJed Brown ierr = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr); 254c4762a1bSJed Brown ctx.Pm = 1.1; 255c4762a1bSJed Brown ierr = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr); 256c4762a1bSJed Brown ctx.tf = 0.1; 257c4762a1bSJed Brown ctx.tcl = 0.2; 258c4762a1bSJed Brown ierr = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr); 259c4762a1bSJed Brown ierr = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr); 260c4762a1bSJed Brown ierr = PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);CHKERRQ(ierr); 261c4762a1bSJed Brown if (ensemble) { 262c4762a1bSJed Brown ctx.tf = -1; 263c4762a1bSJed Brown ctx.tcl = -1; 264c4762a1bSJed Brown } 265c4762a1bSJed Brown 266c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 267c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax); 268c4762a1bSJed Brown u[1] = 1.0; 269c4762a1bSJed Brown ierr = PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);CHKERRQ(ierr); 270c4762a1bSJed Brown n = 2; 271c4762a1bSJed Brown ierr = PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);CHKERRQ(ierr); 272c4762a1bSJed Brown u[0] += du[0]; 273c4762a1bSJed Brown u[1] += du[1]; 274c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 275c4762a1bSJed Brown if (flg1 || flg2) { 276c4762a1bSJed Brown ctx.tf = -1; 277c4762a1bSJed Brown ctx.tcl = -1; 278c4762a1bSJed Brown } 279c4762a1bSJed Brown } 280c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 281c4762a1bSJed Brown 282c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 283c4762a1bSJed Brown Create timestepping solver context 284c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 285c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 286c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 287c4762a1bSJed Brown ierr = TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 288c4762a1bSJed Brown ierr = TSSetType(ts,TSRK);CHKERRQ(ierr); 289c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx);CHKERRQ(ierr); 290c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx);CHKERRQ(ierr); 291c4762a1bSJed Brown ierr = TSCreateQuadratureTS(ts,PETSC_TRUE,&quadts);CHKERRQ(ierr); 292c4762a1bSJed Brown ierr = TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx);CHKERRQ(ierr); 293c4762a1bSJed Brown ierr = TSSetRHSJacobian(quadts,DRDU,DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx);CHKERRQ(ierr); 294c4762a1bSJed Brown ierr = TSSetRHSJacobianP(quadts,DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,&ctx);CHKERRQ(ierr); 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 297c4762a1bSJed Brown Set initial conditions 298c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 299c4762a1bSJed Brown ierr = TSSetSolution(ts,U);CHKERRQ(ierr); 300c4762a1bSJed Brown 301c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 302c4762a1bSJed Brown Save trajectory of solution so that TSAdjointSolve() may be used 303c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 304c4762a1bSJed Brown ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 305c4762a1bSJed Brown 306c4762a1bSJed Brown ierr = MatCreateVecs(A,&lambda[0],NULL);CHKERRQ(ierr); 307c4762a1bSJed Brown /* Set initial conditions for the adjoint integration */ 308c4762a1bSJed Brown ierr = VecGetArray(lambda[0],&y_ptr);CHKERRQ(ierr); 309c4762a1bSJed Brown y_ptr[0] = 0.0; y_ptr[1] = 0.0; 310c4762a1bSJed Brown ierr = VecRestoreArray(lambda[0],&y_ptr);CHKERRQ(ierr); 311c4762a1bSJed Brown 312c4762a1bSJed Brown ierr = MatCreateVecs(Jacp,&mu[0],NULL);CHKERRQ(ierr); 313c4762a1bSJed Brown ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr); 314c4762a1bSJed Brown x_ptr[0] = -1.0; 315c4762a1bSJed Brown ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr); 316c4762a1bSJed Brown ierr = TSSetCostGradients(ts,1,lambda,mu);CHKERRQ(ierr); 317c4762a1bSJed Brown 318c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 319c4762a1bSJed Brown Set solver options 320c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 321c4762a1bSJed Brown ierr = TSSetMaxTime(ts,10.0);CHKERRQ(ierr); 322c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 323c4762a1bSJed Brown ierr = TSSetTimeStep(ts,.01);CHKERRQ(ierr); 324c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 325c4762a1bSJed Brown 326c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 327c4762a1bSJed Brown Solve nonlinear system 328c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 329c4762a1bSJed Brown if (ensemble) { 330c4762a1bSJed Brown for (du[1] = -2.5; du[1] <= .01; du[1] += .1) { 331c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 332c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax); 333c4762a1bSJed Brown u[1] = ctx.omega_s; 334c4762a1bSJed Brown u[0] += du[0]; 335c4762a1bSJed Brown u[1] += du[1]; 336c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 337c4762a1bSJed Brown ierr = TSSetTimeStep(ts,.01);CHKERRQ(ierr); 338c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 339c4762a1bSJed Brown } 340c4762a1bSJed Brown } else { 341c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 342c4762a1bSJed Brown } 343c4762a1bSJed Brown ierr = VecView(U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 344c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 345c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 346c4762a1bSJed Brown 347c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 348c4762a1bSJed Brown Adjoint model starts here 349c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 350c4762a1bSJed Brown /* Set initial conditions for the adjoint integration */ 351c4762a1bSJed Brown ierr = VecGetArray(lambda[0],&y_ptr);CHKERRQ(ierr); 352c4762a1bSJed Brown y_ptr[0] = 0.0; y_ptr[1] = 0.0; 353c4762a1bSJed Brown ierr = VecRestoreArray(lambda[0],&y_ptr);CHKERRQ(ierr); 354c4762a1bSJed Brown 355c4762a1bSJed Brown ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr); 356c4762a1bSJed Brown x_ptr[0] = -1.0; 357c4762a1bSJed Brown ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr); 358c4762a1bSJed Brown 359c4762a1bSJed Brown /* Set RHS JacobianP */ 360c4762a1bSJed Brown ierr = TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&ctx);CHKERRQ(ierr); 361c4762a1bSJed Brown 362c4762a1bSJed Brown ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 363c4762a1bSJed Brown 364c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0] d[Psi(tf)]/d[omega0]\n");CHKERRQ(ierr); 365c4762a1bSJed Brown ierr = VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 366c4762a1bSJed Brown ierr = VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 367c4762a1bSJed Brown ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr); 368c4762a1bSJed Brown ierr = VecView(q,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 369c4762a1bSJed Brown ierr = VecGetArray(q,&x_ptr);CHKERRQ(ierr); 370c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(x_ptr[0]-ctx.Pm));CHKERRQ(ierr); 371c4762a1bSJed Brown ierr = VecRestoreArray(q,&x_ptr);CHKERRQ(ierr); 372c4762a1bSJed Brown 373c4762a1bSJed Brown ierr = ComputeSensiP(lambda[0],mu[0],&ctx);CHKERRQ(ierr); 374c4762a1bSJed Brown 375c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 376c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 377c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 378c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 379c4762a1bSJed Brown ierr = MatDestroy(&Jacp);CHKERRQ(ierr); 380c4762a1bSJed Brown ierr = MatDestroy(&DRDU);CHKERRQ(ierr); 381c4762a1bSJed Brown ierr = MatDestroy(&DRDP);CHKERRQ(ierr); 382c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 383c4762a1bSJed Brown ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr); 384c4762a1bSJed Brown ierr = VecDestroy(&mu[0]);CHKERRQ(ierr); 385c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 386c4762a1bSJed Brown ierr = PetscFinalize(); 387c4762a1bSJed Brown return ierr; 388c4762a1bSJed Brown } 389c4762a1bSJed Brown 390c4762a1bSJed Brown 391c4762a1bSJed Brown /*TEST 392c4762a1bSJed Brown 393c4762a1bSJed Brown build: 394c4762a1bSJed Brown requires: !complex 395c4762a1bSJed Brown 396c4762a1bSJed Brown test: 397c4762a1bSJed Brown args: -viewer_binary_skip_info -ts_adapt_type none 398c4762a1bSJed Brown 399c4762a1bSJed Brown TEST*/ 400