xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 17145e753d664083a23ef2da54d41eaf7d3a7877)
1 /*
2   Code for timestepping with implicit Theta method
3 */
4 #include <petsc-private/tsimpl.h>                /*I   "petscts.h"   I*/
5 #include <petscsnes.h>
6 #include <petscdm.h>
7 #include <petscmat.h>
8 
9 typedef struct {
10   Vec          X,Xdot;                   /* Storage for one stage */
11   Vec          X0;                       /* work vector to store X0 */
12   Vec          affine;                   /* Affine vector needed for residual at beginning of step */
13   Vec          *VecDeltaLam;             /* Increment of the adjoint sensitivity w.r.t IC at stage*/
14   Vec          *VecDeltaMu;              /* Increment of the adjoint sensitivity w.r.t P at stage*/
15   Vec          *VecSensiTemp;            /* Vector to be timed with Jacobian transpose*/
16   PetscBool    extrapolate;
17   PetscBool    endpoint;
18   PetscReal    Theta;
19   PetscReal    stage_time;
20   TSStepStatus status;
21   char         *name;
22   PetscInt     order;
23   PetscReal    ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
24   PetscBool    adapt;  /* use time-step adaptivity ? */
25 } TS_Theta;
26 
27 #undef __FUNCT__
28 #define __FUNCT__ "TSThetaGetX0AndXdot"
29 static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
30 {
31   TS_Theta       *th = (TS_Theta*)ts->data;
32   PetscErrorCode ierr;
33 
34   PetscFunctionBegin;
35   if (X0) {
36     if (dm && dm != ts->dm) {
37       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
38     } else *X0 = ts->vec_sol;
39   }
40   if (Xdot) {
41     if (dm && dm != ts->dm) {
42       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
43     } else *Xdot = th->Xdot;
44   }
45   PetscFunctionReturn(0);
46 }
47 
48 
49 #undef __FUNCT__
50 #define __FUNCT__ "TSThetaRestoreX0AndXdot"
51 static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
52 {
53   PetscErrorCode ierr;
54 
55   PetscFunctionBegin;
56   if (X0) {
57     if (dm && dm != ts->dm) {
58       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
59     }
60   }
61   if (Xdot) {
62     if (dm && dm != ts->dm) {
63       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
64     }
65   }
66   PetscFunctionReturn(0);
67 }
68 
69 #undef __FUNCT__
70 #define __FUNCT__ "DMCoarsenHook_TSTheta"
71 static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
72 {
73 
74   PetscFunctionBegin;
75   PetscFunctionReturn(0);
76 }
77 
78 #undef __FUNCT__
79 #define __FUNCT__ "DMRestrictHook_TSTheta"
80 static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
81 {
82   TS             ts = (TS)ctx;
83   PetscErrorCode ierr;
84   Vec            X0,Xdot,X0_c,Xdot_c;
85 
86   PetscFunctionBegin;
87   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
88   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
89   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
90   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
91   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
92   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
93   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
94   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
95   PetscFunctionReturn(0);
96 }
97 
98 #undef __FUNCT__
99 #define __FUNCT__ "DMSubDomainHook_TSTheta"
100 static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
101 {
102 
103   PetscFunctionBegin;
104   PetscFunctionReturn(0);
105 }
106 
107 #undef __FUNCT__
108 #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta"
109 static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
110 {
111   TS             ts = (TS)ctx;
112   PetscErrorCode ierr;
113   Vec            X0,Xdot,X0_sub,Xdot_sub;
114 
115   PetscFunctionBegin;
116   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
117   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
118 
119   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
120   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
121 
122   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
123   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
124 
125   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
126   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
127   PetscFunctionReturn(0);
128 }
129 
130 #undef __FUNCT__
131 #define __FUNCT__ "TSEvaluateStep_Theta"
132 static PetscErrorCode TSEvaluateStep_Theta(TS ts,PetscInt order,Vec U,PetscBool *done)
133 {
134   PetscErrorCode ierr;
135   TS_Theta       *th = (TS_Theta*)ts->data;
136 
137   PetscFunctionBegin;
138   if (order == 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"No time-step adaptivity implemented for 1st order theta method; Run with -ts_adapt_type none");
139   if (order == th->order) {
140     if (th->endpoint) {
141       ierr = VecCopy(th->X,U);CHKERRQ(ierr);
142     } else {
143       PetscReal shift = 1./(th->Theta*ts->time_step);
144       ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,U,th->X);CHKERRQ(ierr);
145       ierr = VecAXPY(U,ts->time_step,th->Xdot);CHKERRQ(ierr);
146     }
147   } else if (order == th->order-1 && order) {
148     ierr = VecWAXPY(U,ts->time_step,th->Xdot,th->X0);CHKERRQ(ierr);
149   }
150   PetscFunctionReturn(0);
151 }
152 
153 #undef __FUNCT__
154 #define __FUNCT__ "TSRollBack_Theta"
155 static PetscErrorCode TSRollBack_Theta(TS ts)
156 {
157   TS_Theta       *th = (TS_Theta*)ts->data;
158   PetscErrorCode ierr;
159 
160   PetscFunctionBegin;
161   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
162   th->status    = TS_STEP_INCOMPLETE;
163   PetscFunctionReturn(0);
164 }
165 
166 #undef __FUNCT__
167 #define __FUNCT__ "TSStep_Theta"
168 static PetscErrorCode TSStep_Theta(TS ts)
169 {
170   TS_Theta       *th = (TS_Theta*)ts->data;
171   PetscInt       its,lits,reject,next_scheme;
172   PetscReal      next_time_step;
173   TSAdapt        adapt;
174   PetscBool      stageok,accept = PETSC_TRUE;
175   PetscErrorCode ierr;
176 
177   PetscFunctionBegin;
178   th->status = TS_STEP_INCOMPLETE;
179   ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
180   for (reject=0; !ts->reason && th->status != TS_STEP_COMPLETE; ts->reject++) {
181     PetscReal shift = 1./(th->Theta*ts->time_step);
182     th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
183     ierr = TSPreStep(ts);CHKERRQ(ierr);
184     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
185 
186     if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
187       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
188       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
189       ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
190       ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
191     }
192     if (th->extrapolate) {
193       ierr = VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
194     } else {
195       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
196     }
197     ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
198     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
199     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
200     ts->snes_its += its; ts->ksp_its += lits;
201     ierr = TSPostStage(ts,th->stage_time,0,&(th->X));CHKERRQ(ierr);
202     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
203     ierr = TSAdaptCheckStage(adapt,ts,&stageok);CHKERRQ(ierr);
204     if (!stageok) {accept = PETSC_FALSE; goto reject_step;}
205 
206     ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr);
207     th->status = TS_STEP_PENDING;
208     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
209     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
210     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
211     ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr);
212     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
213     if (!accept) {           /* Roll back the current step */
214       ts->ptime += next_time_step; /* This will be undone in rollback */
215       th->status = TS_STEP_INCOMPLETE;
216       ierr = TSRollBack(ts);CHKERRQ(ierr);
217       goto reject_step;
218     }
219 
220     /* This in only used in the context of TSAdjoint. Evaluate cost integral if there is any */
221     if (ts->vec_costintegral) {
222       /* Evolve ts->vec_costintegral to compute integrals */
223       if (th->endpoint) {
224         ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
225         ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(1.-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr);
226       }
227       ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
228       if (th->endpoint) {
229         ierr = VecAXPY(ts->vec_costintegral,ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
230       }else {
231         ierr = VecAXPY(ts->vec_costintegral,ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
232       }
233     }
234 
235     /* ignore next_scheme for now */
236     ts->ptime    += ts->time_step;
237     ts->time_step = next_time_step;
238     ts->steps++;
239     th->status = TS_STEP_COMPLETE;
240     break;
241 
242 reject_step:
243     if (!ts->reason && ++reject > ts->max_reject && ts->max_reject >= 0) {
244       ts->reason = TS_DIVERGED_STEP_REJECTED;
245       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr);
246     }
247     continue;
248   }
249   PetscFunctionReturn(0);
250 }
251 
252 #undef __FUNCT__
253 #define __FUNCT__ "TSAdjointStep_Theta"
254 static PetscErrorCode TSAdjointStep_Theta(TS ts)
255 {
256   TS_Theta            *th = (TS_Theta*)ts->data;
257   Vec                 *VecDeltaLam = th->VecDeltaLam,*VecDeltaMu = th->VecDeltaMu,*VecSensiTemp = th->VecSensiTemp;
258   PetscInt            nadj;
259   PetscErrorCode      ierr;
260   Mat                 J,Jp;
261   KSP                 ksp;
262   PetscReal           shift;
263 
264   PetscFunctionBegin;
265 
266   th->status = TS_STEP_INCOMPLETE;
267   ierr = SNESGetKSP(ts->snes,&ksp);
268   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
269   th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/
270 
271   ierr = TSPreStep(ts);CHKERRQ(ierr);
272 
273   /* Build RHS */
274   if (ts->vec_costintegral) { /* Cost function has an integral  term */
275     if (th->endpoint) {
276       ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr);
277     }else {
278       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
279     }
280   }
281   for (nadj=0; nadj<ts->numcost; nadj++) {
282     ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
283     ierr = VecScale(VecSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr);
284     if (ts->vec_costintegral) {
285       ierr = VecAXPY(VecSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
286     }
287   }
288 
289   /* Build LHS */
290   shift = -1./(th->Theta*ts->time_step);
291   if (th->endpoint) {
292     ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
293   }else {
294     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
295   }
296   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
297 
298   /* Solve LHS X = RHS */
299   for (nadj=0; nadj<ts->numcost; nadj++) {
300     ierr = KSPSolveTranspose(ksp,VecSensiTemp[nadj],VecDeltaLam[nadj]);CHKERRQ(ierr);
301   }
302 
303   /* Update sensitivities, and evaluate integrals if there is any */
304   if(th->endpoint && th->Theta!=1.) { /* two-stage case */
305     shift = -1./((th->Theta-1.)*ts->time_step);
306     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
307     if (ts->vec_costintegral) {
308       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
309     }
310     if (!ts->costintegraldone) {
311       /* Evolve ts->vec_costintegral to compute integrals */
312       ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
313       ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
314       ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
315       ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1.),ts->vec_costintegrand);CHKERRQ(ierr);
316     }
317     for (nadj=0; nadj<ts->numcost; nadj++) {
318       ierr = MatMultTranspose(J,VecDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
319       if (ts->vec_costintegral) {
320         ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
321       }
322       ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr);
323     }
324 
325     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
326       ierr = TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr);
327       for (nadj=0; nadj<ts->numcost; nadj++) {
328         ierr = MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);CHKERRQ(ierr);
329         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecDeltaMu[nadj]);CHKERRQ(ierr);
330       }
331       ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr);
332       for (nadj=0; nadj<ts->numcost; nadj++) {
333         ierr = MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);CHKERRQ(ierr);
334         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecDeltaMu[nadj]);CHKERRQ(ierr);
335       }
336       if (ts->vec_costintegral) {
337         ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr);
338         for (nadj=0; nadj<ts->numcost; nadj++) {
339           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
340         }
341         ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
342         for (nadj=0; nadj<ts->numcost; nadj++) {
343           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr);
344         }
345       }
346     }
347   }else { /* one-stage case */
348     shift = 0.0;
349     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
350     if (ts->vec_costintegral) {
351       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
352     }
353     if (!ts->costintegraldone) {
354       /* Evolve ts->vec_costintegral to compute integrals */
355       ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
356       ierr = VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
357     }
358      /* When th->endpoint is true and th->Theta==1 (beuler method), the Jacobian is supposed to be evaluated at ts->ptime like this:
359     if(th->endpoint) {
360       ierr  = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
361     }
362     but ts->ptime and ts->vec_sol have the same values as th->stage_time and th->X in this case. So the code is simplified here.
363     */
364     for (nadj=0; nadj<ts->numcost; nadj++) {
365       ierr = MatMultTranspose(J,VecDeltaLam[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
366       ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecSensiTemp[nadj]);CHKERRQ(ierr);
367       if (ts->vec_costintegral) {
368         ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
369       }
370     }
371     if (ts->vecs_sensip) {
372       ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr);
373       for (nadj=0; nadj<ts->numcost; nadj++) {
374         ierr = MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);CHKERRQ(ierr);
375         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecDeltaMu[nadj]);CHKERRQ(ierr);
376       }
377       if (ts->vec_costintegral) {
378         ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
379         for (nadj=0; nadj<ts->numcost; nadj++) {
380           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
381         }
382       }
383     }
384   }
385 
386   ts->ptime += ts->time_step;
387   ts->steps++;
388   th->status = TS_STEP_COMPLETE;
389   PetscFunctionReturn(0);
390 }
391 
392 #undef __FUNCT__
393 #define __FUNCT__ "TSInterpolate_Theta"
394 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
395 {
396   TS_Theta       *th   = (TS_Theta*)ts->data;
397   PetscReal      alpha = t - ts->ptime;
398   PetscErrorCode ierr;
399 
400   PetscFunctionBegin;
401   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
402   if (th->endpoint) alpha *= th->Theta;
403   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
404   PetscFunctionReturn(0);
405 }
406 
407 /*------------------------------------------------------------*/
408 #undef __FUNCT__
409 #define __FUNCT__ "TSReset_Theta"
410 static PetscErrorCode TSReset_Theta(TS ts)
411 {
412   TS_Theta       *th = (TS_Theta*)ts->data;
413   PetscErrorCode ierr;
414 
415   PetscFunctionBegin;
416   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
417   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
418   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
419   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
420   ierr = VecDestroyVecs(ts->numcost,&th->VecDeltaLam);CHKERRQ(ierr);
421   ierr = VecDestroyVecs(ts->numcost,&th->VecDeltaMu);CHKERRQ(ierr);
422   ierr = VecDestroyVecs(ts->numcost,&th->VecSensiTemp);CHKERRQ(ierr);
423   PetscFunctionReturn(0);
424 }
425 
426 #undef __FUNCT__
427 #define __FUNCT__ "TSDestroy_Theta"
428 static PetscErrorCode TSDestroy_Theta(TS ts)
429 {
430   PetscErrorCode ierr;
431 
432   PetscFunctionBegin;
433   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
434   ierr = PetscFree(ts->data);CHKERRQ(ierr);
435   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
436   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
437   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
438   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
439   PetscFunctionReturn(0);
440 }
441 
442 /*
443   This defines the nonlinear equation that is to be solved with SNES
444   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
445 */
446 #undef __FUNCT__
447 #define __FUNCT__ "SNESTSFormFunction_Theta"
448 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
449 {
450   TS_Theta       *th = (TS_Theta*)ts->data;
451   PetscErrorCode ierr;
452   Vec            X0,Xdot;
453   DM             dm,dmsave;
454   PetscReal      shift = 1./(th->Theta*ts->time_step);
455 
456   PetscFunctionBegin;
457   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
458   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
459   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
460   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
461 
462   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
463   dmsave = ts->dm;
464   ts->dm = dm;
465   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
466   ts->dm = dmsave;
467   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
468   PetscFunctionReturn(0);
469 }
470 
471 #undef __FUNCT__
472 #define __FUNCT__ "SNESTSFormJacobian_Theta"
473 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
474 {
475   TS_Theta       *th = (TS_Theta*)ts->data;
476   PetscErrorCode ierr;
477   Vec            Xdot;
478   DM             dm,dmsave;
479   PetscReal      shift = 1./(th->Theta*ts->time_step);
480 
481   PetscFunctionBegin;
482   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
483 
484   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
485   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
486 
487   dmsave = ts->dm;
488   ts->dm = dm;
489   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
490   ts->dm = dmsave;
491   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
492   PetscFunctionReturn(0);
493 }
494 
495 #undef __FUNCT__
496 #define __FUNCT__ "TSSetUp_Theta"
497 static PetscErrorCode TSSetUp_Theta(TS ts)
498 {
499   TS_Theta       *th = (TS_Theta*)ts->data;
500   PetscErrorCode ierr;
501   SNES           snes;
502   TSAdapt        adapt;
503   DM             dm;
504 
505   PetscFunctionBegin;
506   if (!th->X) {
507     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
508   }
509   if (!th->Xdot) {
510     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
511   }
512   if (!th->X0) {
513     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
514   }
515   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
516   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
517   if (dm) {
518     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
519     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
520   }
521   if (th->Theta == 0.5 && th->endpoint) th->order = 2;
522   else th->order = 1;
523 
524   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
525   if (!th->adapt) {
526     ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr);
527   }
528   PetscFunctionReturn(0);
529 }
530 /*------------------------------------------------------------*/
531 
532 #undef __FUNCT__
533 #define __FUNCT__ "TSAdjointSetUp_Theta"
534 static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
535 {
536   TS_Theta       *th = (TS_Theta*)ts->data;
537   PetscErrorCode ierr;
538 
539   PetscFunctionBegin;
540   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecDeltaLam);CHKERRQ(ierr);
541   if(ts->vecs_sensip) {
542     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecDeltaMu);CHKERRQ(ierr);
543   }
544   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecSensiTemp);CHKERRQ(ierr);
545   PetscFunctionReturn(0);
546 }
547 /*------------------------------------------------------------*/
548 
549 #undef __FUNCT__
550 #define __FUNCT__ "TSSetFromOptions_Theta"
551 static PetscErrorCode TSSetFromOptions_Theta(PetscOptions *PetscOptionsObject,TS ts)
552 {
553   TS_Theta       *th = (TS_Theta*)ts->data;
554   PetscErrorCode ierr;
555 
556   PetscFunctionBegin;
557   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
558   {
559     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
560     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
561     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr);
562     ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr);
563     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
564   }
565   ierr = PetscOptionsTail();CHKERRQ(ierr);
566   PetscFunctionReturn(0);
567 }
568 
569 #undef __FUNCT__
570 #define __FUNCT__ "TSView_Theta"
571 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
572 {
573   TS_Theta       *th = (TS_Theta*)ts->data;
574   PetscBool      iascii;
575   PetscErrorCode ierr;
576 
577   PetscFunctionBegin;
578   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
579   if (iascii) {
580     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
581     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
582   }
583   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
584   PetscFunctionReturn(0);
585 }
586 
587 #undef __FUNCT__
588 #define __FUNCT__ "TSThetaGetTheta_Theta"
589 PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
590 {
591   TS_Theta *th = (TS_Theta*)ts->data;
592 
593   PetscFunctionBegin;
594   *theta = th->Theta;
595   PetscFunctionReturn(0);
596 }
597 
598 #undef __FUNCT__
599 #define __FUNCT__ "TSThetaSetTheta_Theta"
600 PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
601 {
602   TS_Theta *th = (TS_Theta*)ts->data;
603 
604   PetscFunctionBegin;
605   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
606   th->Theta = theta;
607   PetscFunctionReturn(0);
608 }
609 
610 #undef __FUNCT__
611 #define __FUNCT__ "TSThetaGetEndpoint_Theta"
612 PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
613 {
614   TS_Theta *th = (TS_Theta*)ts->data;
615 
616   PetscFunctionBegin;
617   *endpoint = th->endpoint;
618   PetscFunctionReturn(0);
619 }
620 
621 #undef __FUNCT__
622 #define __FUNCT__ "TSThetaSetEndpoint_Theta"
623 PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
624 {
625   TS_Theta *th = (TS_Theta*)ts->data;
626 
627   PetscFunctionBegin;
628   th->endpoint = flg;
629   PetscFunctionReturn(0);
630 }
631 
632 #if defined(PETSC_HAVE_COMPLEX)
633 #undef __FUNCT__
634 #define __FUNCT__ "TSComputeLinearStability_Theta"
635 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
636 {
637   PetscComplex z   = xr + xi*PETSC_i,f;
638   TS_Theta     *th = (TS_Theta*)ts->data;
639   const PetscReal one = 1.0;
640 
641   PetscFunctionBegin;
642   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
643   *yr = PetscRealPartComplex(f);
644   *yi = PetscImaginaryPartComplex(f);
645   PetscFunctionReturn(0);
646 }
647 #endif
648 
649 #undef __FUNCT__
650 #define __FUNCT__ "TSGetStages_Theta"
651 static PetscErrorCode  TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
652 {
653   TS_Theta     *th = (TS_Theta*)ts->data;
654 
655   PetscFunctionBegin;
656   *ns = 1;
657   if(Y) {
658     *Y  = &(th->X);
659   }
660   PetscFunctionReturn(0);
661 }
662 
663 /* ------------------------------------------------------------ */
664 /*MC
665       TSTHETA - DAE solver using the implicit Theta method
666 
667    Level: beginner
668 
669    Options Database:
670       -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
671       -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable)
672       -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
673 
674    Notes:
675 $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
676 $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
677 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
678 
679 
680 
681    This method can be applied to DAE.
682 
683    This method is cast as a 1-stage implicit Runge-Kutta method.
684 
685 .vb
686   Theta | Theta
687   -------------
688         |  1
689 .ve
690 
691    For the default Theta=0.5, this is also known as the implicit midpoint rule.
692 
693    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
694 
695 .vb
696   0 | 0         0
697   1 | 1-Theta   Theta
698   -------------------
699     | 1-Theta   Theta
700 .ve
701 
702    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
703 
704    To apply a diagonally implicit RK method to DAE, the stage formula
705 
706 $  Y_i = X + h sum_j a_ij Y'_j
707 
708    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
709 
710 .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
711 
712 M*/
713 #undef __FUNCT__
714 #define __FUNCT__ "TSCreate_Theta"
715 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
716 {
717   TS_Theta       *th;
718   PetscErrorCode ierr;
719 
720   PetscFunctionBegin;
721   ts->ops->reset           = TSReset_Theta;
722   ts->ops->destroy         = TSDestroy_Theta;
723   ts->ops->view            = TSView_Theta;
724   ts->ops->setup           = TSSetUp_Theta;
725   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
726   ts->ops->step            = TSStep_Theta;
727   ts->ops->interpolate     = TSInterpolate_Theta;
728   ts->ops->evaluatestep    = TSEvaluateStep_Theta;
729   ts->ops->rollback        = TSRollBack_Theta;
730   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
731   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
732   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
733 #if defined(PETSC_HAVE_COMPLEX)
734   ts->ops->linearstability = TSComputeLinearStability_Theta;
735 #endif
736   ts->ops->getstages       = TSGetStages_Theta;
737   ts->ops->adjointstep     = TSAdjointStep_Theta;
738 
739   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
740   ts->data = (void*)th;
741 
742   th->extrapolate = PETSC_FALSE;
743   th->Theta       = 0.5;
744   th->ccfl        = 1.0;
745   th->adapt       = PETSC_FALSE;
746   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
747   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
748   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
749   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
750   PetscFunctionReturn(0);
751 }
752 
753 #undef __FUNCT__
754 #define __FUNCT__ "TSThetaGetTheta"
755 /*@
756   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
757 
758   Not Collective
759 
760   Input Parameter:
761 .  ts - timestepping context
762 
763   Output Parameter:
764 .  theta - stage abscissa
765 
766   Note:
767   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
768 
769   Level: Advanced
770 
771 .seealso: TSThetaSetTheta()
772 @*/
773 PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
774 {
775   PetscErrorCode ierr;
776 
777   PetscFunctionBegin;
778   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
779   PetscValidPointer(theta,2);
780   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
781   PetscFunctionReturn(0);
782 }
783 
784 #undef __FUNCT__
785 #define __FUNCT__ "TSThetaSetTheta"
786 /*@
787   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
788 
789   Not Collective
790 
791   Input Parameter:
792 +  ts - timestepping context
793 -  theta - stage abscissa
794 
795   Options Database:
796 .  -ts_theta_theta <theta>
797 
798   Level: Intermediate
799 
800 .seealso: TSThetaGetTheta()
801 @*/
802 PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
803 {
804   PetscErrorCode ierr;
805 
806   PetscFunctionBegin;
807   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
808   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
809   PetscFunctionReturn(0);
810 }
811 
812 #undef __FUNCT__
813 #define __FUNCT__ "TSThetaGetEndpoint"
814 /*@
815   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
816 
817   Not Collective
818 
819   Input Parameter:
820 .  ts - timestepping context
821 
822   Output Parameter:
823 .  endpoint - PETSC_TRUE when using the endpoint variant
824 
825   Level: Advanced
826 
827 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
828 @*/
829 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
830 {
831   PetscErrorCode ierr;
832 
833   PetscFunctionBegin;
834   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
835   PetscValidPointer(endpoint,2);
836   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
837   PetscFunctionReturn(0);
838 }
839 
840 #undef __FUNCT__
841 #define __FUNCT__ "TSThetaSetEndpoint"
842 /*@
843   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
844 
845   Not Collective
846 
847   Input Parameter:
848 +  ts - timestepping context
849 -  flg - PETSC_TRUE to use the endpoint variant
850 
851   Options Database:
852 .  -ts_theta_endpoint <flg>
853 
854   Level: Intermediate
855 
856 .seealso: TSTHETA, TSCN
857 @*/
858 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
859 {
860   PetscErrorCode ierr;
861 
862   PetscFunctionBegin;
863   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
864   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
865   PetscFunctionReturn(0);
866 }
867 
868 /*
869  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
870  * The creation functions for these specializations are below.
871  */
872 
873 #undef __FUNCT__
874 #define __FUNCT__ "TSView_BEuler"
875 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
876 {
877   PetscErrorCode ierr;
878 
879   PetscFunctionBegin;
880   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
881   PetscFunctionReturn(0);
882 }
883 
884 /*MC
885       TSBEULER - ODE solver using the implicit backward Euler method
886 
887   Level: beginner
888 
889   Notes:
890   TSBEULER is equivalent to TSTHETA with Theta=1.0
891 
892 $  -ts_type theta -ts_theta_theta 1.
893 
894 .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
895 
896 M*/
897 #undef __FUNCT__
898 #define __FUNCT__ "TSCreate_BEuler"
899 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
900 {
901   PetscErrorCode ierr;
902 
903   PetscFunctionBegin;
904   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
905   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
906   ts->ops->view = TSView_BEuler;
907   PetscFunctionReturn(0);
908 }
909 
910 #undef __FUNCT__
911 #define __FUNCT__ "TSView_CN"
912 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
913 {
914   PetscErrorCode ierr;
915 
916   PetscFunctionBegin;
917   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
918   PetscFunctionReturn(0);
919 }
920 
921 /*MC
922       TSCN - ODE solver using the implicit Crank-Nicolson method.
923 
924   Level: beginner
925 
926   Notes:
927   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
928 
929 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
930 
931 .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
932 
933 M*/
934 #undef __FUNCT__
935 #define __FUNCT__ "TSCreate_CN"
936 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
937 {
938   PetscErrorCode ierr;
939 
940   PetscFunctionBegin;
941   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
942   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
943   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
944   ts->ops->view = TSView_CN;
945   PetscFunctionReturn(0);
946 }
947