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