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