1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Adjoint and tangent linear sensitivity analysis of the basic equation for generator stability analysis.\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown /*F 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown \begin{eqnarray} 7*c4762a1bSJed Brown \frac{d \theta}{dt} = \omega_b (\omega - \omega_s) 8*c4762a1bSJed Brown \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\ 9*c4762a1bSJed Brown \end{eqnarray} 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown F*/ 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown /* 14*c4762a1bSJed Brown This code demonstrate the sensitivity analysis interface to a system of ordinary differential equations with discontinuities. 15*c4762a1bSJed Brown It computes the sensitivities of an integral cost function 16*c4762a1bSJed Brown \int c*max(0,\theta(t)-u_s)^beta dt 17*c4762a1bSJed Brown w.r.t. initial conditions and the parameter P_m. 18*c4762a1bSJed Brown Backward Euler method is used for time integration. 19*c4762a1bSJed Brown The discontinuities are detected with TSEvent. 20*c4762a1bSJed Brown */ 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown #include <petscts.h> 23*c4762a1bSJed Brown #include "ex3.h" 24*c4762a1bSJed Brown 25*c4762a1bSJed Brown int main(int argc,char **argv) 26*c4762a1bSJed Brown { 27*c4762a1bSJed Brown TS ts,quadts; /* ODE integrator */ 28*c4762a1bSJed Brown Vec U; /* solution will be stored here */ 29*c4762a1bSJed Brown PetscErrorCode ierr; 30*c4762a1bSJed Brown PetscMPIInt size; 31*c4762a1bSJed Brown PetscInt n = 2; 32*c4762a1bSJed Brown AppCtx ctx; 33*c4762a1bSJed Brown PetscScalar *u; 34*c4762a1bSJed Brown PetscReal du[2] = {0.0,0.0}; 35*c4762a1bSJed Brown PetscBool ensemble = PETSC_FALSE,flg1,flg2; 36*c4762a1bSJed Brown PetscReal ftime; 37*c4762a1bSJed Brown PetscInt steps; 38*c4762a1bSJed Brown PetscScalar *x_ptr,*y_ptr,*s_ptr; 39*c4762a1bSJed Brown Vec lambda[1],q,mu[1]; 40*c4762a1bSJed Brown PetscInt direction[2]; 41*c4762a1bSJed Brown PetscBool terminate[2]; 42*c4762a1bSJed Brown Mat qgrad; 43*c4762a1bSJed Brown Mat sp; /* Forward sensitivity matrix */ 44*c4762a1bSJed Brown SAMethod sa; 45*c4762a1bSJed Brown 46*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 47*c4762a1bSJed Brown Initialize program 48*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 49*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 50*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 51*c4762a1bSJed Brown if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs"); 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 54*c4762a1bSJed Brown Create necessary matrix and vectors 55*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 56*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&ctx.Jac);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = MatSetSizes(ctx.Jac,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = MatSetType(ctx.Jac,MATDENSE);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = MatSetFromOptions(ctx.Jac);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = MatSetUp(ctx.Jac);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = MatCreateVecs(ctx.Jac,&U,NULL);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&ctx.Jacp);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = MatSetSizes(ctx.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = MatSetFromOptions(ctx.Jacp);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = MatSetUp(ctx.Jacp);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&ctx.DRDP);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = MatSetUp(ctx.DRDP);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&ctx.DRDU);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = MatSetUp(ctx.DRDU);CHKERRQ(ierr); 70*c4762a1bSJed Brown 71*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 72*c4762a1bSJed Brown Set runtime options 73*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 74*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr); 75*c4762a1bSJed Brown { 76*c4762a1bSJed Brown ctx.beta = 2; 77*c4762a1bSJed Brown ctx.c = 10000.0; 78*c4762a1bSJed Brown ctx.u_s = 1.0; 79*c4762a1bSJed Brown ctx.omega_s = 1.0; 80*c4762a1bSJed Brown ctx.omega_b = 120.0*PETSC_PI; 81*c4762a1bSJed Brown ctx.H = 5.0; 82*c4762a1bSJed Brown ierr = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr); 83*c4762a1bSJed Brown ctx.D = 5.0; 84*c4762a1bSJed Brown ierr = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr); 85*c4762a1bSJed Brown ctx.E = 1.1378; 86*c4762a1bSJed Brown ctx.V = 1.0; 87*c4762a1bSJed Brown ctx.X = 0.545; 88*c4762a1bSJed Brown ctx.Pmax = ctx.E*ctx.V/ctx.X; 89*c4762a1bSJed Brown ctx.Pmax_ini = ctx.Pmax; 90*c4762a1bSJed Brown ierr = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr); 91*c4762a1bSJed Brown ctx.Pm = 1.1; 92*c4762a1bSJed Brown ierr = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr); 93*c4762a1bSJed Brown ctx.tf = 0.1; 94*c4762a1bSJed Brown ctx.tcl = 0.2; 95*c4762a1bSJed Brown ierr = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);CHKERRQ(ierr); 98*c4762a1bSJed Brown if (ensemble) { 99*c4762a1bSJed Brown ctx.tf = -1; 100*c4762a1bSJed Brown ctx.tcl = -1; 101*c4762a1bSJed Brown } 102*c4762a1bSJed Brown 103*c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 104*c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax); 105*c4762a1bSJed Brown u[1] = 1.0; 106*c4762a1bSJed Brown ierr = PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);CHKERRQ(ierr); 107*c4762a1bSJed Brown n = 2; 108*c4762a1bSJed Brown ierr = PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);CHKERRQ(ierr); 109*c4762a1bSJed Brown u[0] += du[0]; 110*c4762a1bSJed Brown u[1] += du[1]; 111*c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 112*c4762a1bSJed Brown if (flg1 || flg2) { 113*c4762a1bSJed Brown ctx.tf = -1; 114*c4762a1bSJed Brown ctx.tcl = -1; 115*c4762a1bSJed Brown } 116*c4762a1bSJed Brown sa = SA_ADJ; 117*c4762a1bSJed Brown ierr = PetscOptionsEnum("-sa_method","Sensitivity analysis method (adj or tlm)","",SAMethods,(PetscEnum)sa,(PetscEnum*)&sa,NULL);CHKERRQ(ierr); 118*c4762a1bSJed Brown } 119*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 120*c4762a1bSJed Brown 121*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 122*c4762a1bSJed Brown Create timestepping solver context 123*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 124*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr); 127*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx);CHKERRQ(ierr); 128*c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,ctx.Jac,ctx.Jac,(TSRHSJacobian)RHSJacobian,&ctx);CHKERRQ(ierr); 129*c4762a1bSJed Brown 130*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 131*c4762a1bSJed Brown Set initial conditions 132*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 133*c4762a1bSJed Brown ierr = TSSetSolution(ts,U);CHKERRQ(ierr); 134*c4762a1bSJed Brown 135*c4762a1bSJed Brown /* Set RHS JacobianP */ 136*c4762a1bSJed Brown ierr = TSSetRHSJacobianP(ts,ctx.Jacp,RHSJacobianP,&ctx);CHKERRQ(ierr); 137*c4762a1bSJed Brown 138*c4762a1bSJed Brown ierr = TSCreateQuadratureTS(ts,PETSC_FALSE,&quadts);CHKERRQ(ierr); 139*c4762a1bSJed Brown ierr = TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx);CHKERRQ(ierr); 140*c4762a1bSJed Brown ierr = TSSetRHSJacobian(quadts,ctx.DRDU,ctx.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx);CHKERRQ(ierr); 141*c4762a1bSJed Brown ierr = TSSetRHSJacobianP(quadts,ctx.DRDP,DRDPJacobianTranspose,&ctx);CHKERRQ(ierr); 142*c4762a1bSJed Brown if (sa == SA_ADJ) { 143*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 144*c4762a1bSJed Brown Save trajectory of solution so that TSAdjointSolve() may be used 145*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 146*c4762a1bSJed Brown ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 147*c4762a1bSJed Brown ierr = MatCreateVecs(ctx.Jac,&lambda[0],NULL);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = MatCreateVecs(ctx.Jacp,&mu[0],NULL);CHKERRQ(ierr); 149*c4762a1bSJed Brown ierr = TSSetCostGradients(ts,1,lambda,mu);CHKERRQ(ierr); 150*c4762a1bSJed Brown } 151*c4762a1bSJed Brown 152*c4762a1bSJed Brown if (sa == SA_TLM) { 153*c4762a1bSJed Brown PetscScalar val[2]; 154*c4762a1bSJed Brown PetscInt row[]={0,1},col[]={0}; 155*c4762a1bSJed Brown 156*c4762a1bSJed Brown ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&qgrad);CHKERRQ(ierr); 157*c4762a1bSJed Brown ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&sp);CHKERRQ(ierr); 158*c4762a1bSJed Brown ierr = TSForwardSetSensitivities(ts,1,sp);CHKERRQ(ierr); 159*c4762a1bSJed Brown ierr = TSForwardSetSensitivities(quadts,1,qgrad);CHKERRQ(ierr); 160*c4762a1bSJed Brown val[0] = 1./PetscSqrtScalar(1.-(ctx.Pm/ctx.Pmax)*(ctx.Pm/ctx.Pmax))/ctx.Pmax; 161*c4762a1bSJed Brown val[1] = 0.0; 162*c4762a1bSJed Brown ierr = MatSetValues(sp,2,row,1,col,val,INSERT_VALUES);CHKERRQ(ierr); 163*c4762a1bSJed Brown ierr = MatAssemblyBegin(sp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = MatAssemblyEnd(sp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 165*c4762a1bSJed Brown } 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 168*c4762a1bSJed Brown Set solver options 169*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 170*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,1.0);CHKERRQ(ierr); 171*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 172*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.03125);CHKERRQ(ierr); 173*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 174*c4762a1bSJed Brown 175*c4762a1bSJed Brown direction[0] = direction[1] = 1; 176*c4762a1bSJed Brown terminate[0] = terminate[1] = PETSC_FALSE; 177*c4762a1bSJed Brown 178*c4762a1bSJed Brown ierr = TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)&ctx);CHKERRQ(ierr); 179*c4762a1bSJed Brown 180*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 181*c4762a1bSJed Brown Solve nonlinear system 182*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 183*c4762a1bSJed Brown if (ensemble) { 184*c4762a1bSJed Brown for (du[1] = -2.5; du[1] <= .01; du[1] += .1) { 185*c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 186*c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax); 187*c4762a1bSJed Brown u[1] = ctx.omega_s; 188*c4762a1bSJed Brown u[0] += du[0]; 189*c4762a1bSJed Brown u[1] += du[1]; 190*c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 191*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.03125);CHKERRQ(ierr); 192*c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 193*c4762a1bSJed Brown } 194*c4762a1bSJed Brown } else { 195*c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 196*c4762a1bSJed Brown } 197*c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 198*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 199*c4762a1bSJed Brown 200*c4762a1bSJed Brown if (sa == SA_ADJ) { 201*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 202*c4762a1bSJed Brown Adjoint model starts here 203*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 204*c4762a1bSJed Brown /* Set initial conditions for the adjoint integration */ 205*c4762a1bSJed Brown ierr = VecGetArray(lambda[0],&y_ptr);CHKERRQ(ierr); 206*c4762a1bSJed Brown y_ptr[0] = 0.0; y_ptr[1] = 0.0; 207*c4762a1bSJed Brown ierr = VecRestoreArray(lambda[0],&y_ptr);CHKERRQ(ierr); 208*c4762a1bSJed Brown 209*c4762a1bSJed Brown ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr); 210*c4762a1bSJed Brown x_ptr[0] = 0.0; 211*c4762a1bSJed Brown ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr); 212*c4762a1bSJed Brown 213*c4762a1bSJed Brown ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 214*c4762a1bSJed Brown 215*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n lambda: d[Psi(tf)]/d[phi0] d[Psi(tf)]/d[omega0]\n");CHKERRQ(ierr); 216*c4762a1bSJed Brown ierr = VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 217*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n mu: d[Psi(tf)]/d[pm]\n");CHKERRQ(ierr); 218*c4762a1bSJed Brown ierr = VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 219*c4762a1bSJed Brown ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr); 220*c4762a1bSJed Brown ierr = VecGetArray(q,&x_ptr);CHKERRQ(ierr); 221*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(x_ptr[0]-ctx.Pm));CHKERRQ(ierr); 222*c4762a1bSJed Brown ierr = VecRestoreArray(q,&x_ptr);CHKERRQ(ierr); 223*c4762a1bSJed Brown ierr = ComputeSensiP(lambda[0],mu[0],&ctx);CHKERRQ(ierr); 224*c4762a1bSJed Brown ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr); 225*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n gradient=%g\n",(double)x_ptr[0]);CHKERRQ(ierr); 226*c4762a1bSJed Brown ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr); 227*c4762a1bSJed Brown ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr); 228*c4762a1bSJed Brown ierr = VecDestroy(&mu[0]);CHKERRQ(ierr); 229*c4762a1bSJed Brown } 230*c4762a1bSJed Brown if (sa == SA_TLM) { 231*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n trajectory sensitivity: d[phi(tf)]/d[pm] d[omega(tf)]/d[pm]\n");CHKERRQ(ierr); 232*c4762a1bSJed Brown ierr = MatView(sp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 233*c4762a1bSJed Brown ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr); 234*c4762a1bSJed Brown ierr = VecGetArray(q,&s_ptr);CHKERRQ(ierr); 235*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(s_ptr[0]-ctx.Pm));CHKERRQ(ierr); 236*c4762a1bSJed Brown ierr = VecRestoreArray(q,&s_ptr);CHKERRQ(ierr); 237*c4762a1bSJed Brown ierr = MatDenseGetArray(qgrad,&s_ptr);CHKERRQ(ierr); 238*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n gradient=%g\n",(double)s_ptr[0]);CHKERRQ(ierr); 239*c4762a1bSJed Brown ierr = MatDenseRestoreArray(qgrad,&s_ptr);CHKERRQ(ierr); 240*c4762a1bSJed Brown ierr = MatDestroy(&qgrad);CHKERRQ(ierr); 241*c4762a1bSJed Brown ierr = MatDestroy(&sp);CHKERRQ(ierr); 242*c4762a1bSJed Brown } 243*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 244*c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 245*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 246*c4762a1bSJed Brown ierr = MatDestroy(&ctx.Jac);CHKERRQ(ierr); 247*c4762a1bSJed Brown ierr = MatDestroy(&ctx.Jacp);CHKERRQ(ierr); 248*c4762a1bSJed Brown ierr = MatDestroy(&ctx.DRDU);CHKERRQ(ierr); 249*c4762a1bSJed Brown ierr = MatDestroy(&ctx.DRDP);CHKERRQ(ierr); 250*c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 251*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 252*c4762a1bSJed Brown ierr = PetscFinalize(); 253*c4762a1bSJed Brown return ierr; 254*c4762a1bSJed Brown } 255*c4762a1bSJed Brown 256*c4762a1bSJed Brown 257*c4762a1bSJed Brown /*TEST 258*c4762a1bSJed Brown 259*c4762a1bSJed Brown build: 260*c4762a1bSJed Brown requires: !complex !single 261*c4762a1bSJed Brown 262*c4762a1bSJed Brown test: 263*c4762a1bSJed Brown args: -sa_method adj -viewer_binary_skip_info -ts_type cn -pc_type lu 264*c4762a1bSJed Brown 265*c4762a1bSJed Brown test: 266*c4762a1bSJed Brown suffix: 2 267*c4762a1bSJed Brown args: -sa_method tlm -ts_type cn -pc_type lu 268*c4762a1bSJed Brown 269*c4762a1bSJed Brown test: 270*c4762a1bSJed Brown suffix: 3 271*c4762a1bSJed Brown args: -sa_method adj -ts_type rk -ts_rk_type 2a -ts_adapt_type dsp 272*c4762a1bSJed Brown 273*c4762a1bSJed Brown TEST*/ 274