xref: /petsc/src/ts/tutorials/power_grid/ex3sa.c (revision 3c633725528e547aaaa9b672a746f5d686a276e1)
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;
50ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
51*3c633725SBarry 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     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&ctx.Jac);CHKERRQ(ierr);
57c4762a1bSJed Brown   ierr = MatSetSizes(ctx.Jac,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
58c4762a1bSJed Brown   ierr = MatSetType(ctx.Jac,MATDENSE);CHKERRQ(ierr);
59c4762a1bSJed Brown   ierr = MatSetFromOptions(ctx.Jac);CHKERRQ(ierr);
60c4762a1bSJed Brown   ierr = MatSetUp(ctx.Jac);CHKERRQ(ierr);
61c4762a1bSJed Brown   ierr = MatCreateVecs(ctx.Jac,&U,NULL);CHKERRQ(ierr);
62c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&ctx.Jacp);CHKERRQ(ierr);
63c4762a1bSJed Brown   ierr = MatSetSizes(ctx.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr);
64c4762a1bSJed Brown   ierr = MatSetFromOptions(ctx.Jacp);CHKERRQ(ierr);
65c4762a1bSJed Brown   ierr = MatSetUp(ctx.Jacp);CHKERRQ(ierr);
66c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&ctx.DRDP);CHKERRQ(ierr);
67c4762a1bSJed Brown   ierr = MatSetUp(ctx.DRDP);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&ctx.DRDU);CHKERRQ(ierr);
69c4762a1bSJed Brown   ierr = MatSetUp(ctx.DRDU);CHKERRQ(ierr);
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;
82c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr);
83c4762a1bSJed Brown     ctx.D       = 5.0;
84c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr);
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;
90c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr);
91c4762a1bSJed Brown     ctx.Pm      = 1.1;
92c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr);
93c4762a1bSJed Brown     ctx.tf      = 0.1;
94c4762a1bSJed Brown     ctx.tcl     = 0.2;
95c4762a1bSJed Brown     ierr        = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr);
96c4762a1bSJed Brown     ierr        = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr);
97c4762a1bSJed Brown     ierr        = PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);CHKERRQ(ierr);
98c4762a1bSJed Brown     if (ensemble) {
99c4762a1bSJed Brown       ctx.tf      = -1;
100c4762a1bSJed Brown       ctx.tcl     = -1;
101c4762a1bSJed Brown     }
102c4762a1bSJed Brown 
103c4762a1bSJed Brown     ierr = VecGetArray(U,&u);CHKERRQ(ierr);
104c4762a1bSJed Brown     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
105c4762a1bSJed Brown     u[1] = 1.0;
106c4762a1bSJed Brown     ierr = PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);CHKERRQ(ierr);
107c4762a1bSJed Brown     n    = 2;
108c4762a1bSJed Brown     ierr = PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);CHKERRQ(ierr);
109c4762a1bSJed Brown     u[0] += du[0];
110c4762a1bSJed Brown     u[1] += du[1];
111c4762a1bSJed Brown     ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
112c4762a1bSJed Brown     if (flg1 || flg2) {
113c4762a1bSJed Brown       ctx.tf      = -1;
114c4762a1bSJed Brown       ctx.tcl     = -1;
115c4762a1bSJed Brown     }
116c4762a1bSJed Brown     sa = SA_ADJ;
117c4762a1bSJed Brown     ierr = PetscOptionsEnum("-sa_method","Sensitivity analysis method (adj or tlm)","",SAMethods,(PetscEnum)sa,(PetscEnum*)&sa,NULL);CHKERRQ(ierr);
118c4762a1bSJed Brown   }
119c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122c4762a1bSJed Brown      Create timestepping solver context
123c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
125c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
127c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx);CHKERRQ(ierr);
128c4762a1bSJed Brown   ierr = TSSetRHSJacobian(ts,ctx.Jac,ctx.Jac,(TSRHSJacobian)RHSJacobian,&ctx);CHKERRQ(ierr);
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131c4762a1bSJed Brown      Set initial conditions
132c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133c4762a1bSJed Brown   ierr = TSSetSolution(ts,U);CHKERRQ(ierr);
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /*   Set RHS JacobianP */
136c4762a1bSJed Brown   ierr = TSSetRHSJacobianP(ts,ctx.Jacp,RHSJacobianP,&ctx);CHKERRQ(ierr);
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   ierr = TSCreateQuadratureTS(ts,PETSC_FALSE,&quadts);CHKERRQ(ierr);
139c4762a1bSJed Brown   ierr = TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx);CHKERRQ(ierr);
140c4762a1bSJed Brown   ierr = TSSetRHSJacobian(quadts,ctx.DRDU,ctx.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx);CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr = TSSetRHSJacobianP(quadts,ctx.DRDP,DRDPJacobianTranspose,&ctx);CHKERRQ(ierr);
142c4762a1bSJed Brown   if (sa == SA_ADJ) {
143c4762a1bSJed Brown     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144c4762a1bSJed Brown       Save trajectory of solution so that TSAdjointSolve() may be used
145c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146c4762a1bSJed Brown     ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
147c4762a1bSJed Brown     ierr = MatCreateVecs(ctx.Jac,&lambda[0],NULL);CHKERRQ(ierr);
148c4762a1bSJed Brown     ierr = MatCreateVecs(ctx.Jacp,&mu[0],NULL);CHKERRQ(ierr);
149c4762a1bSJed Brown     ierr = TSSetCostGradients(ts,1,lambda,mu);CHKERRQ(ierr);
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 
156c4762a1bSJed Brown     ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&qgrad);CHKERRQ(ierr);
157c4762a1bSJed Brown     ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&sp);CHKERRQ(ierr);
158c4762a1bSJed Brown     ierr = TSForwardSetSensitivities(ts,1,sp);CHKERRQ(ierr);
159c4762a1bSJed Brown     ierr = TSForwardSetSensitivities(quadts,1,qgrad);CHKERRQ(ierr);
160c4762a1bSJed Brown     val[0] = 1./PetscSqrtScalar(1.-(ctx.Pm/ctx.Pmax)*(ctx.Pm/ctx.Pmax))/ctx.Pmax;
161c4762a1bSJed Brown     val[1] = 0.0;
162c4762a1bSJed Brown     ierr = MatSetValues(sp,2,row,1,col,val,INSERT_VALUES);CHKERRQ(ierr);
163c4762a1bSJed Brown     ierr = MatAssemblyBegin(sp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
164c4762a1bSJed Brown     ierr = MatAssemblyEnd(sp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
165c4762a1bSJed Brown   }
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168c4762a1bSJed Brown      Set solver options
169c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,1.0);CHKERRQ(ierr);
171c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
172c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,0.03125);CHKERRQ(ierr);
173c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   direction[0] = direction[1] = 1;
176c4762a1bSJed Brown   terminate[0] = terminate[1] = PETSC_FALSE;
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   ierr = TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)&ctx);CHKERRQ(ierr);
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) {
185c4762a1bSJed Brown       ierr = VecGetArray(U,&u);CHKERRQ(ierr);
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];
190c4762a1bSJed Brown       ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
191c4762a1bSJed Brown       ierr = TSSetTimeStep(ts,0.03125);CHKERRQ(ierr);
192c4762a1bSJed Brown       ierr = TSSolve(ts,U);CHKERRQ(ierr);
193c4762a1bSJed Brown     }
194c4762a1bSJed Brown   } else {
195c4762a1bSJed Brown     ierr = TSSolve(ts,U);CHKERRQ(ierr);
196c4762a1bSJed Brown   }
197c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
198c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
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 */
205c4762a1bSJed Brown     ierr = VecGetArray(lambda[0],&y_ptr);CHKERRQ(ierr);
206c4762a1bSJed Brown     y_ptr[0] = 0.0; y_ptr[1] = 0.0;
207c4762a1bSJed Brown     ierr = VecRestoreArray(lambda[0],&y_ptr);CHKERRQ(ierr);
208c4762a1bSJed Brown 
209c4762a1bSJed Brown     ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr);
210c4762a1bSJed Brown     x_ptr[0] = 0.0;
211c4762a1bSJed Brown     ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr);
212c4762a1bSJed Brown 
213c4762a1bSJed Brown     ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
214c4762a1bSJed Brown 
215c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n lambda: d[Psi(tf)]/d[phi0]  d[Psi(tf)]/d[omega0]\n");CHKERRQ(ierr);
216c4762a1bSJed Brown     ierr = VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
217c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n mu: d[Psi(tf)]/d[pm]\n");CHKERRQ(ierr);
218c4762a1bSJed Brown     ierr = VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
219c4762a1bSJed Brown     ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr);
220c4762a1bSJed Brown     ierr = VecGetArray(q,&x_ptr);CHKERRQ(ierr);
221c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(x_ptr[0]-ctx.Pm));CHKERRQ(ierr);
222c4762a1bSJed Brown     ierr = VecRestoreArray(q,&x_ptr);CHKERRQ(ierr);
223c4762a1bSJed Brown     ierr = ComputeSensiP(lambda[0],mu[0],&ctx);CHKERRQ(ierr);
224c4762a1bSJed Brown     ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr);
225c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n gradient=%g\n",(double)x_ptr[0]);CHKERRQ(ierr);
226c4762a1bSJed Brown     ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr);
227c4762a1bSJed Brown     ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr);
228c4762a1bSJed Brown     ierr = VecDestroy(&mu[0]);CHKERRQ(ierr);
229c4762a1bSJed Brown   }
230c4762a1bSJed Brown   if (sa == SA_TLM) {
231c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n trajectory sensitivity: d[phi(tf)]/d[pm]  d[omega(tf)]/d[pm]\n");CHKERRQ(ierr);
232c4762a1bSJed Brown     ierr = MatView(sp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
233c4762a1bSJed Brown     ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr);
234c4762a1bSJed Brown     ierr = VecGetArray(q,&s_ptr);CHKERRQ(ierr);
235c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(s_ptr[0]-ctx.Pm));CHKERRQ(ierr);
236c4762a1bSJed Brown     ierr = VecRestoreArray(q,&s_ptr);CHKERRQ(ierr);
237c4762a1bSJed Brown     ierr = MatDenseGetArray(qgrad,&s_ptr);CHKERRQ(ierr);
238c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n gradient=%g\n",(double)s_ptr[0]);CHKERRQ(ierr);
239c4762a1bSJed Brown     ierr = MatDenseRestoreArray(qgrad,&s_ptr);CHKERRQ(ierr);
240c4762a1bSJed Brown     ierr = MatDestroy(&qgrad);CHKERRQ(ierr);
241c4762a1bSJed Brown     ierr = MatDestroy(&sp);CHKERRQ(ierr);
242c4762a1bSJed Brown   }
243c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
245c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
246c4762a1bSJed Brown   ierr = MatDestroy(&ctx.Jac);CHKERRQ(ierr);
247c4762a1bSJed Brown   ierr = MatDestroy(&ctx.Jacp);CHKERRQ(ierr);
248c4762a1bSJed Brown   ierr = MatDestroy(&ctx.DRDU);CHKERRQ(ierr);
249c4762a1bSJed Brown   ierr = MatDestroy(&ctx.DRDP);CHKERRQ(ierr);
250c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
251c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
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