Actual source code: ts.c

petsc-3.3-p7 2013-05-11
  2: #include <petsc-private/tsimpl.h>        /*I "petscts.h"  I*/
  3: #include <petscdmshell.h>

  5: /* Logging support */
  6: PetscClassId  TS_CLASSID;
  7: PetscLogEvent  TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;

 11: /*
 12:   TSSetTypeFromOptions - Sets the type of ts from user options.

 14:   Collective on TS

 16:   Input Parameter:
 17: . ts - The ts

 19:   Level: intermediate

 21: .keywords: TS, set, options, database, type
 22: .seealso: TSSetFromOptions(), TSSetType()
 23: */
 24: static PetscErrorCode TSSetTypeFromOptions(TS ts)
 25: {
 26:   PetscBool      opt;
 27:   const char     *defaultType;
 28:   char           typeName[256];

 32:   if (((PetscObject)ts)->type_name) {
 33:     defaultType = ((PetscObject)ts)->type_name;
 34:   } else {
 35:     defaultType = TSEULER;
 36:   }

 38:   if (!TSRegisterAllCalled) {TSRegisterAll(PETSC_NULL);}
 39:   PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);
 40:   if (opt) {
 41:     TSSetType(ts, typeName);
 42:   } else {
 43:     TSSetType(ts, defaultType);
 44:   }
 45:   return(0);
 46: }

 50: /*@
 51:    TSSetFromOptions - Sets various TS parameters from user options.

 53:    Collective on TS

 55:    Input Parameter:
 56: .  ts - the TS context obtained from TSCreate()

 58:    Options Database Keys:
 59: +  -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP
 60: .  -ts_max_steps maxsteps - maximum number of time-steps to take
 61: .  -ts_final_time time - maximum time to compute to
 62: .  -ts_dt dt - initial time step
 63: .  -ts_monitor - print information at each timestep
 64: -  -ts_monitor_draw - plot information at each timestep

 66:    Level: beginner

 68: .keywords: TS, timestep, set, options, database

 70: .seealso: TSGetType()
 71: @*/
 72: PetscErrorCode  TSSetFromOptions(TS ts)
 73: {
 74:   PetscBool      opt,flg;
 76:   PetscViewer    monviewer;
 77:   char           monfilename[PETSC_MAX_PATH_LEN];
 78:   SNES           snes;
 79:   TSAdapt        adapt;

 83:   PetscObjectOptionsBegin((PetscObject)ts);
 84:     /* Handle TS type options */
 85:     TSSetTypeFromOptions(ts);

 87:     /* Handle generic TS options */
 88:     PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);
 89:     PetscOptionsReal("-ts_final_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);
 90:     PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,PETSC_NULL);
 91:     PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&ts->time_step,PETSC_NULL);
 92:     opt = ts->exact_final_time == PETSC_DECIDE ? PETSC_FALSE : (PetscBool)ts->exact_final_time;
 93:     PetscOptionsBool("-ts_exact_final_time","Interpolate output to stop exactly at the final time","TSSetExactFinalTime",opt,&opt,&flg);
 94:     if (flg) {TSSetExactFinalTime(ts,opt);}
 95:     PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","TSSetMaxSNESFailures",ts->max_snes_failures,&ts->max_snes_failures,PETSC_NULL);
 96:     PetscOptionsInt("-ts_max_reject","Maximum number of step rejections before step fails","TSSetMaxStepRejections",ts->max_reject,&ts->max_reject,PETSC_NULL);
 97:     PetscOptionsBool("-ts_error_if_step_fails","Error if no step succeeds","TSSetErrorIfStepFails",ts->errorifstepfailed,&ts->errorifstepfailed,PETSC_NULL);
 98:     PetscOptionsReal("-ts_rtol","Relative tolerance for local truncation error","TSSetTolerances",ts->rtol,&ts->rtol,PETSC_NULL);
 99:     PetscOptionsReal("-ts_atol","Absolute tolerance for local truncation error","TSSetTolerances",ts->atol,&ts->atol,PETSC_NULL);

101:     /* Monitor options */
102:     PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
103:     if (flg) {
104:       PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);
105:       TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
106:     }
107:     PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
108:     if (flg) {PetscPythonMonitorSet((PetscObject)ts,monfilename);}

110:     opt  = PETSC_FALSE;
111:     PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);
112:     if (opt) {
113:       TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);
114:     }
115:     opt  = PETSC_FALSE;
116:     PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);
117:     if (opt) {
118:       void *ctx;
119:       TSMonitorSolutionCreate(ts,PETSC_NULL,&ctx);
120:       TSMonitorSet(ts,TSMonitorSolution,ctx,TSMonitorSolutionDestroy);
121:     }
122:     opt  = PETSC_FALSE;
123:     PetscOptionsString("-ts_monitor_solution_binary","Save each solution to a binary file","TSMonitorSolutionBinary",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
124:     if (flg) {
125:       PetscViewer ctx;
126:       if (monfilename[0]) {
127:         PetscViewerBinaryOpen(((PetscObject)ts)->comm,monfilename,FILE_MODE_WRITE,&ctx);
128:       } else {
129:         ctx = PETSC_VIEWER_BINARY_(((PetscObject)ts)->comm);
130:       }
131:       TSMonitorSet(ts,TSMonitorSolutionBinary,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
132:     }
133:     opt  = PETSC_FALSE;
134:     PetscOptionsString("-ts_monitor_solution_vtk","Save each time step to a binary file, use filename-%%03D.vts","TSMonitorSolutionVTK",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
135:     if (flg) {
136:       const char *ptr,*ptr2;
137:       char *filetemplate;
138:       if (!monfilename[0]) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
139:       /* Do some cursory validation of the input. */
140:       PetscStrstr(monfilename,"%",(char**)&ptr);
141:       if (!ptr) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
142:       for (ptr++ ; ptr && *ptr; ptr++) {
143:         PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);
144:         if (!ptr2 && (*ptr < '0' || '9' < *ptr)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Invalid file template argument to -ts_monitor_solution_vtk, should look like filename-%%03D.vts");
145:         if (ptr2) break;
146:       }
147:       PetscStrallocpy(monfilename,&filetemplate);
148:       TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);
149:     }

151:     TSGetAdapt(ts,&adapt);
152:     TSAdaptSetFromOptions(adapt);

154:     TSGetSNES(ts,&snes);
155:     if (ts->problem_type == TS_LINEAR) {SNESSetType(snes,SNESKSPONLY);}

157:     /* Handle specific TS options */
158:     if (ts->ops->setfromoptions) {
159:       (*ts->ops->setfromoptions)(ts);
160:     }

162:     /* process any options handlers added with PetscObjectAddOptionsHandler() */
163:     PetscObjectProcessOptionsHandlers((PetscObject)ts);
164:   PetscOptionsEnd();
165:   return(0);
166: }

171: /*@
172:    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
173:       set with TSSetRHSJacobian().

175:    Collective on TS and Vec

177:    Input Parameters:
178: +  ts - the TS context
179: .  t - current timestep
180: -  x - input vector

182:    Output Parameters:
183: +  A - Jacobian matrix
184: .  B - optional preconditioning matrix
185: -  flag - flag indicating matrix structure

187:    Notes:
188:    Most users should not need to explicitly call this routine, as it
189:    is used internally within the nonlinear solvers.

191:    See KSPSetOperators() for important information about setting the
192:    flag parameter.

194:    Level: developer

196: .keywords: SNES, compute, Jacobian, matrix

198: .seealso:  TSSetRHSJacobian(), KSPSetOperators()
199: @*/
200: PetscErrorCode  TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
201: {
203:   PetscInt Xstate;

209:   PetscObjectStateQuery((PetscObject)X,&Xstate);
210:   if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == X && ts->rhsjacobian.Xstate == Xstate))) {
211:     *flg = ts->rhsjacobian.mstructure;
212:     return(0);
213:   }

215:   if (!ts->userops->rhsjacobian && !ts->userops->ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");

217:   if (ts->userops->rhsjacobian) {
218:     PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);
219:     *flg = DIFFERENT_NONZERO_PATTERN;
220:     PetscStackPush("TS user Jacobian function");
221:     (*ts->userops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);
222:     PetscStackPop;
223:     PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);
224:     /* make sure user returned a correct Jacobian and preconditioner */
227:   } else {
228:     MatZeroEntries(*A);
229:     if (*A != *B) {MatZeroEntries(*B);}
230:     *flg = SAME_NONZERO_PATTERN;
231:   }
232:   ts->rhsjacobian.time = t;
233:   ts->rhsjacobian.X = X;
234:   PetscObjectStateQuery((PetscObject)X,&ts->rhsjacobian.Xstate);
235:   ts->rhsjacobian.mstructure = *flg;
236:   return(0);
237: }

241: /*@
242:    TSComputeRHSFunction - Evaluates the right-hand-side function. 

244:    Collective on TS and Vec

246:    Input Parameters:
247: +  ts - the TS context
248: .  t - current time
249: -  x - state vector

251:    Output Parameter:
252: .  y - right hand side

254:    Note:
255:    Most users should not need to explicitly call this routine, as it
256:    is used internally within the nonlinear solvers.

258:    Level: developer

260: .keywords: TS, compute

262: .seealso: TSSetRHSFunction(), TSComputeIFunction()
263: @*/
264: PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
265: {


273:   if (!ts->userops->rhsfunction && !ts->userops->ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");

275:   PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);
276:   if (ts->userops->rhsfunction) {
277:     PetscStackPush("TS user right-hand-side function");
278:     (*ts->userops->rhsfunction)(ts,t,x,y,ts->funP);
279:     PetscStackPop;
280:   } else {
281:     VecZeroEntries(y);
282:   }

284:   PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);
285:   return(0);
286: }

290: static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
291: {
292:   Vec            F;

296:   *Frhs = PETSC_NULL;
297:   TSGetIFunction(ts,&F,PETSC_NULL,PETSC_NULL);
298:   if (!ts->Frhs) {
299:     VecDuplicate(F,&ts->Frhs);
300:   }
301:   *Frhs = ts->Frhs;
302:   return(0);
303: }

307: static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
308: {
309:   Mat            A,B;

313:   TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);
314:   if (Arhs) {
315:     if (!ts->Arhs) {
316:       MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);
317:     }
318:     *Arhs = ts->Arhs;
319:   }
320:   if (Brhs) {
321:     if (!ts->Brhs) {
322:       MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);
323:     }
324:     *Brhs = ts->Brhs;
325:   }
326:   return(0);
327: }

331: /*@
332:    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0

334:    Collective on TS and Vec

336:    Input Parameters:
337: +  ts - the TS context
338: .  t - current time
339: .  X - state vector
340: .  Xdot - time derivative of state vector
341: -  imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate

343:    Output Parameter:
344: .  Y - right hand side

346:    Note:
347:    Most users should not need to explicitly call this routine, as it
348:    is used internally within the nonlinear solvers.

350:    If the user did did not write their equations in implicit form, this
351:    function recasts them in implicit form.

353:    Level: developer

355: .keywords: TS, compute

357: .seealso: TSSetIFunction(), TSComputeRHSFunction()
358: @*/
359: PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex)
360: {


369:   if (!ts->userops->rhsfunction && !ts->userops->ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");

371:   PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);
372:   if (ts->userops->ifunction) {
373:     PetscStackPush("TS user implicit function");
374:     (*ts->userops->ifunction)(ts,t,X,Xdot,Y,ts->funP);
375:     PetscStackPop;
376:   }
377:   if (imex) {
378:     if (!ts->userops->ifunction) {
379:       VecCopy(Xdot,Y);
380:     }
381:   } else if (ts->userops->rhsfunction) {
382:     if (ts->userops->ifunction) {
383:       Vec Frhs;
384:       TSGetRHSVec_Private(ts,&Frhs);
385:       TSComputeRHSFunction(ts,t,X,Frhs);
386:       VecAXPY(Y,-1,Frhs);
387:     } else {
388:       TSComputeRHSFunction(ts,t,X,Y);
389:       VecAYPX(Y,-1,Xdot);
390:     }
391:   }
392:   PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);
393:   return(0);
394: }

398: /*@
399:    TSComputeIJacobian - Evaluates the Jacobian of the DAE

401:    Collective on TS and Vec

403:    Input
404:       Input Parameters:
405: +  ts - the TS context
406: .  t - current timestep
407: .  X - state vector
408: .  Xdot - time derivative of state vector
409: .  shift - shift to apply, see note below
410: -  imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate

412:    Output Parameters:
413: +  A - Jacobian matrix
414: .  B - optional preconditioning matrix
415: -  flag - flag indicating matrix structure

417:    Notes:
418:    If F(t,X,Xdot)=0 is the DAE, the required Jacobian is

420:    dF/dX + shift*dF/dXdot

422:    Most users should not need to explicitly call this routine, as it
423:    is used internally within the nonlinear solvers.

425:    Level: developer

427: .keywords: TS, compute, Jacobian, matrix

429: .seealso:  TSSetIJacobian()
430: @*/
431: PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex)
432: {
433:   PetscInt Xstate, Xdotstate;

445:   PetscObjectStateQuery((PetscObject)X,&Xstate);
446:   PetscObjectStateQuery((PetscObject)Xdot,&Xdotstate);
447:   if (ts->ijacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->ijacobian.X == X && ts->ijacobian.Xstate == Xstate && ts->ijacobian.Xdot == Xdot && ts->ijacobian.Xdotstate == Xdotstate && ts->ijacobian.imex == imex))) {
448:     *flg = ts->ijacobian.mstructure;
449:     MatScale(*A, shift / ts->ijacobian.shift);
450:     return(0);
451:   }

453:   if (!ts->userops->rhsjacobian && !ts->userops->ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");

455:   *flg = SAME_NONZERO_PATTERN;  /* In case we're solving a linear problem in which case it wouldn't get initialized below. */
456:   PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);
457:   if (ts->userops->ijacobian) {
458:     *flg = DIFFERENT_NONZERO_PATTERN;
459:     PetscStackPush("TS user implicit Jacobian");
460:     (*ts->userops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);
461:     PetscStackPop;
462:     /* make sure user returned a correct Jacobian and preconditioner */
465:   }
466:   if (imex) {
467:     if (!ts->userops->ijacobian) {  /* system was written as Xdot = F(t,X) */
468:       MatZeroEntries(*A);
469:       MatShift(*A,shift);
470:       if (*A != *B) {
471:         MatZeroEntries(*B);
472:         MatShift(*B,shift);
473:       }
474:       *flg = SAME_PRECONDITIONER;
475:     }
476:   } else {
477:     if (!ts->userops->ijacobian) {
478:       TSComputeRHSJacobian(ts,t,X,A,B,flg);
479:       MatScale(*A,-1);
480:       MatShift(*A,shift);
481:       if (*A != *B) {
482:         MatScale(*B,-1);
483:         MatShift(*B,shift);
484:       }
485:     } else if (ts->userops->rhsjacobian) {
486:       Mat Arhs,Brhs;
487:       MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN;
488:       TSGetRHSMats_Private(ts,&Arhs,&Brhs);
489:       TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);
490:       axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN;
491:       MatAXPY(*A,-1,Arhs,axpy);
492:       if (*A != *B) {
493:         MatAXPY(*B,-1,Brhs,axpy);
494:       }
495:       *flg = PetscMin(*flg,flg2);
496:     }
497:   }

499:   ts->ijacobian.time = t;
500:   ts->ijacobian.X = X;
501:   ts->ijacobian.Xdot = Xdot;
502:   PetscObjectStateQuery((PetscObject)X,&ts->ijacobian.Xstate);
503:   PetscObjectStateQuery((PetscObject)Xdot,&ts->ijacobian.Xdotstate);
504:   ts->ijacobian.shift = shift;
505:   ts->ijacobian.imex = imex;
506:   ts->ijacobian.mstructure = *flg;
507:   PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);
508:   return(0);
509: }

513: /*@C
514:     TSSetRHSFunction - Sets the routine for evaluating the function,
515:     F(t,u), where U_t = F(t,u).

517:     Logically Collective on TS

519:     Input Parameters:
520: +   ts - the TS context obtained from TSCreate()
521: .   r - vector to put the computed right hand side (or PETSC_NULL to have it created)
522: .   f - routine for evaluating the right-hand-side function
523: -   ctx - [optional] user-defined context for private data for the 
524:           function evaluation routine (may be PETSC_NULL)

526:     Calling sequence of func:
527: $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);

529: +   t - current timestep
530: .   u - input vector
531: .   F - function vector
532: -   ctx - [optional] user-defined function context 

534:     Level: beginner

536: .keywords: TS, timestep, set, right-hand-side, function

538: .seealso: TSSetRHSJacobian(), TSSetIJacobian()
539: @*/
540: PetscErrorCode  TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
541: {
543:   SNES           snes;
544:   Vec            ralloc = PETSC_NULL;

549:   if (f)   ts->userops->rhsfunction = f;
550:   if (ctx) ts->funP                 = ctx;
551:   TSGetSNES(ts,&snes);
552:   if (!r && !ts->dm && ts->vec_sol) {
553:     VecDuplicate(ts->vec_sol,&ralloc);
554:     r = ralloc;
555:   }
556:   SNESSetFunction(snes,r,SNESTSFormFunction,ts);
557:   VecDestroy(&ralloc);
558:   return(0);
559: }

563: /*@C
564:    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
565:    where U_t = F(U,t), as well as the location to store the matrix.

567:    Logically Collective on TS

569:    Input Parameters:
570: +  ts  - the TS context obtained from TSCreate()
571: .  A   - Jacobian matrix
572: .  B   - preconditioner matrix (usually same as A)
573: .  f   - the Jacobian evaluation routine
574: -  ctx - [optional] user-defined context for private data for the 
575:          Jacobian evaluation routine (may be PETSC_NULL)

577:    Calling sequence of func:
578: $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);

580: +  t - current timestep
581: .  u - input vector
582: .  A - matrix A, where U_t = A(t)u
583: .  B - preconditioner matrix, usually the same as A
584: .  flag - flag indicating information about the preconditioner matrix
585:           structure (same as flag in KSPSetOperators())
586: -  ctx - [optional] user-defined context for matrix evaluation routine

588:    Notes: 
589:    See KSPSetOperators() for important information about setting the flag
590:    output parameter in the routine func().  Be sure to read this information!

592:    The routine func() takes Mat * as the matrix arguments rather than Mat.  
593:    This allows the matrix evaluation routine to replace A and/or B with a 
594:    completely new matrix structure (not just different matrix elements)
595:    when appropriate, for instance, if the nonzero structure is changing
596:    throughout the global iterations.

598:    Level: beginner
599:    
600: .keywords: TS, timestep, set, right-hand-side, Jacobian

602: .seealso: SNESDefaultComputeJacobianColor(), TSSetRHSFunction()

604: @*/
605: PetscErrorCode  TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx)
606: {
608:   SNES           snes;


617:   if (f)   ts->userops->rhsjacobian = f;
618:   if (ctx) ts->jacP                 = ctx;
619:   TSGetSNES(ts,&snes);
620:   if (!ts->userops->ijacobian) {
621:     SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);
622:   }
623:   if (A) {
624:     PetscObjectReference((PetscObject)A);
625:     MatDestroy(&ts->Arhs);
626:     ts->Arhs = A;
627:   }
628:   if (B) {
629:     PetscObjectReference((PetscObject)B);
630:     MatDestroy(&ts->Brhs);
631:     ts->Brhs = B;
632:   }
633:   return(0);
634: }


639: /*@C
640:    TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved.

642:    Logically Collective on TS

644:    Input Parameters:
645: +  ts  - the TS context obtained from TSCreate()
646: .  r   - vector to hold the residual (or PETSC_NULL to have it created internally)
647: .  f   - the function evaluation routine
648: -  ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL)

650:    Calling sequence of f:
651: $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);

653: +  t   - time at step/stage being solved
654: .  u   - state vector
655: .  u_t - time derivative of state vector
656: .  F   - function vector
657: -  ctx - [optional] user-defined context for matrix evaluation routine

659:    Important:
660:    The user MUST call either this routine, TSSetRHSFunction().  This routine must be used when not solving an ODE, for example a DAE.

662:    Level: beginner

664: .keywords: TS, timestep, set, DAE, Jacobian

666: .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian()
667: @*/
668: PetscErrorCode  TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx)
669: {
671:   SNES           snes;
672:   Vec            resalloc = PETSC_NULL;

677:   if (f)   ts->userops->ifunction = f;
678:   if (ctx) ts->funP           = ctx;
679:   TSGetSNES(ts,&snes);
680:   if (!res && !ts->dm && ts->vec_sol) {
681:     VecDuplicate(ts->vec_sol,&resalloc);
682:     res = resalloc;
683:   }
684:   SNESSetFunction(snes,res,SNESTSFormFunction,ts);
685:   VecDestroy(&resalloc);
686:   return(0);
687: }

691: /*@C
692:    TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it.

694:    Not Collective

696:    Input Parameter:
697: .  ts - the TS context

699:    Output Parameter:
700: +  r - vector to hold residual (or PETSC_NULL)
701: .  func - the function to compute residual (or PETSC_NULL)
702: -  ctx - the function context (or PETSC_NULL)

704:    Level: advanced

706: .keywords: TS, nonlinear, get, function

708: .seealso: TSSetIFunction(), SNESGetFunction()
709: @*/
710: PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
711: {
713:   SNES snes;

717:   TSGetSNES(ts,&snes);
718:   SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);
719:   if (func) *func = ts->userops->ifunction;
720:   if (ctx)  *ctx  = ts->funP;
721:   return(0);
722: }

726: /*@C
727:    TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it.

729:    Not Collective

731:    Input Parameter:
732: .  ts - the TS context

734:    Output Parameter:
735: +  r - vector to hold computed right hand side (or PETSC_NULL)
736: .  func - the function to compute right hand side (or PETSC_NULL)
737: -  ctx - the function context (or PETSC_NULL)

739:    Level: advanced

741: .keywords: TS, nonlinear, get, function

743: .seealso: TSSetRhsfunction(), SNESGetFunction()
744: @*/
745: PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
746: {
748:   SNES snes;

752:   TSGetSNES(ts,&snes);
753:   SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);
754:   if (func) *func = ts->userops->rhsfunction;
755:   if (ctx)  *ctx  = ts->funP;
756:   return(0);
757: }

761: /*@C
762:    TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
763:         you provided with TSSetIFunction().  

765:    Logically Collective on TS

767:    Input Parameters:
768: +  ts  - the TS context obtained from TSCreate()
769: .  A   - Jacobian matrix
770: .  B   - preconditioning matrix for A (may be same as A)
771: .  f   - the Jacobian evaluation routine
772: -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL)

774:    Calling sequence of f:
775: $  f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx);

777: +  t    - time at step/stage being solved
778: .  U    - state vector
779: .  U_t  - time derivative of state vector
780: .  a    - shift
781: .  A    - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
782: .  B    - preconditioning matrix for A, may be same as A
783: .  flag - flag indicating information about the preconditioner matrix
784:           structure (same as flag in KSPSetOperators())
785: -  ctx  - [optional] user-defined context for matrix evaluation routine

787:    Notes:
788:    The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve.

790:    The matrix dF/dU + a*dF/dU_t you provide turns out to be 
791:    the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
792:    The time integrator internally approximates U_t by W+a*U where the positive "shift"
793:    a and vector W depend on the integration method, step size, and past states. For example with 
794:    the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
795:    W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt

797:    Level: beginner

799: .keywords: TS, timestep, DAE, Jacobian

801: .seealso: TSSetIFunction(), TSSetRHSJacobian(), SNESDefaultComputeJacobianColor(), SNESDefaultComputeJacobian()

803: @*/
804: PetscErrorCode  TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx)
805: {
807:   SNES           snes;

815:   if (f)   ts->userops->ijacobian = f;
816:   if (ctx) ts->jacP           = ctx;
817:   TSGetSNES(ts,&snes);
818:   SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);
819:   return(0);
820: }

824: /*@C
825:     TSView - Prints the TS data structure.

827:     Collective on TS

829:     Input Parameters:
830: +   ts - the TS context obtained from TSCreate()
831: -   viewer - visualization context

833:     Options Database Key:
834: .   -ts_view - calls TSView() at end of TSStep()

836:     Notes:
837:     The available visualization contexts include
838: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
839: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
840:          output where only the first processor opens
841:          the file.  All other processors send their
842:          data to the first processor to print.

844:     The user can open an alternative visualization context with
845:     PetscViewerASCIIOpen() - output to a specified file.

847:     Level: beginner

849: .keywords: TS, timestep, view

851: .seealso: PetscViewerASCIIOpen()
852: @*/
853: PetscErrorCode  TSView(TS ts,PetscViewer viewer)
854: {
856:   const TSType   type;
857:   PetscBool      iascii,isstring,isundials;

861:   if (!viewer) {
862:     PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);
863:   }

867:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
868:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
869:   if (iascii) {
870:     PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");
871:     PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);
872:     PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);
873:     if (ts->problem_type == TS_NONLINEAR) {
874:       PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->snes_its);
875:       PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solve failures=%D\n",ts->num_snes_failures);
876:     }
877:     PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->ksp_its);
878:     PetscViewerASCIIPrintf(viewer,"  total number of rejected steps=%D\n",ts->reject);
879:     if (ts->ops->view) {
880:       PetscViewerASCIIPushTab(viewer);
881:       (*ts->ops->view)(ts,viewer);
882:       PetscViewerASCIIPopTab(viewer);
883:     }
884:   } else if (isstring) {
885:     TSGetType(ts,&type);
886:     PetscViewerStringSPrintf(viewer," %-7.7s",type);
887:   }
888:   PetscViewerASCIIPushTab(viewer);
889:   PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);
890:   PetscViewerASCIIPopTab(viewer);
891:   return(0);
892: }


897: /*@
898:    TSSetApplicationContext - Sets an optional user-defined context for
899:    the timesteppers.

901:    Logically Collective on TS

903:    Input Parameters:
904: +  ts - the TS context obtained from TSCreate()
905: -  usrP - optional user context

907:    Level: intermediate

909: .keywords: TS, timestep, set, application, context

911: .seealso: TSGetApplicationContext()
912: @*/
913: PetscErrorCode  TSSetApplicationContext(TS ts,void *usrP)
914: {
917:   ts->user = usrP;
918:   return(0);
919: }

923: /*@
924:     TSGetApplicationContext - Gets the user-defined context for the 
925:     timestepper.

927:     Not Collective

929:     Input Parameter:
930: .   ts - the TS context obtained from TSCreate()

932:     Output Parameter:
933: .   usrP - user context

935:     Level: intermediate

937: .keywords: TS, timestep, get, application, context

939: .seealso: TSSetApplicationContext()
940: @*/
941: PetscErrorCode  TSGetApplicationContext(TS ts,void *usrP)
942: {
945:   *(void**)usrP = ts->user;
946:   return(0);
947: }

951: /*@
952:    TSGetTimeStepNumber - Gets the number of time steps completed.

954:    Not Collective

956:    Input Parameter:
957: .  ts - the TS context obtained from TSCreate()

959:    Output Parameter:
960: .  iter - number of steps completed so far

962:    Level: intermediate

964: .keywords: TS, timestep, get, iteration, number
965: .seealso: TSGetTime(), TSGetTimeStep(), TSSetPreStep(), TSSetPreStage(), TSSetPostStep()
966: @*/
967: PetscErrorCode  TSGetTimeStepNumber(TS ts,PetscInt* iter)
968: {
972:   *iter = ts->steps;
973:   return(0);
974: }

978: /*@
979:    TSSetInitialTimeStep - Sets the initial timestep to be used,
980:    as well as the initial time.

982:    Logically Collective on TS

984:    Input Parameters:
985: +  ts - the TS context obtained from TSCreate()
986: .  initial_time - the initial time
987: -  time_step - the size of the timestep

989:    Level: intermediate

991: .seealso: TSSetTimeStep(), TSGetTimeStep()

993: .keywords: TS, set, initial, timestep
994: @*/
995: PetscErrorCode  TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
996: {

1001:   TSSetTimeStep(ts,time_step);
1002:   TSSetTime(ts,initial_time);
1003:   return(0);
1004: }

1008: /*@
1009:    TSSetTimeStep - Allows one to reset the timestep at any time,
1010:    useful for simple pseudo-timestepping codes.

1012:    Logically Collective on TS

1014:    Input Parameters:
1015: +  ts - the TS context obtained from TSCreate()
1016: -  time_step - the size of the timestep

1018:    Level: intermediate

1020: .seealso: TSSetInitialTimeStep(), TSGetTimeStep()

1022: .keywords: TS, set, timestep
1023: @*/
1024: PetscErrorCode  TSSetTimeStep(TS ts,PetscReal time_step)
1025: {
1029:   ts->time_step = time_step;
1030:   return(0);
1031: }

1035: /*@
1036:    TSSetExactFinalTime - Determines whether to interpolate solution to the
1037:       exact final time requested by the user or just returns it at the final time
1038:       it computed.

1040:   Logically Collective on TS

1042:    Input Parameter:
1043: +   ts - the time-step context
1044: -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE

1046:    Level: beginner

1048: .seealso: TSSetDuration()
1049: @*/
1050: PetscErrorCode  TSSetExactFinalTime(TS ts,PetscBool flg)
1051: {

1056:   ts->exact_final_time = flg;
1057:   return(0);
1058: }

1062: /*@
1063:    TSGetTimeStep - Gets the current timestep size.

1065:    Not Collective

1067:    Input Parameter:
1068: .  ts - the TS context obtained from TSCreate()

1070:    Output Parameter:
1071: .  dt - the current timestep size

1073:    Level: intermediate

1075: .seealso: TSSetInitialTimeStep(), TSGetTimeStep()

1077: .keywords: TS, get, timestep
1078: @*/
1079: PetscErrorCode  TSGetTimeStep(TS ts,PetscReal* dt)
1080: {
1084:   *dt = ts->time_step;
1085:   return(0);
1086: }

1090: /*@
1091:    TSGetSolution - Returns the solution at the present timestep. It
1092:    is valid to call this routine inside the function that you are evaluating
1093:    in order to move to the new timestep. This vector not changed until
1094:    the solution at the next timestep has been calculated.

1096:    Not Collective, but Vec returned is parallel if TS is parallel

1098:    Input Parameter:
1099: .  ts - the TS context obtained from TSCreate()

1101:    Output Parameter:
1102: .  v - the vector containing the solution

1104:    Level: intermediate

1106: .seealso: TSGetTimeStep()

1108: .keywords: TS, timestep, get, solution
1109: @*/
1110: PetscErrorCode  TSGetSolution(TS ts,Vec *v)
1111: {
1115:   *v = ts->vec_sol;
1116:   return(0);
1117: }

1119: /* ----- Routines to initialize and destroy a timestepper ---- */
1122: /*@
1123:   TSSetProblemType - Sets the type of problem to be solved.

1125:   Not collective

1127:   Input Parameters:
1128: + ts   - The TS
1129: - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1130: .vb
1131:          U_t = A U    
1132:          U_t = A(t) U 
1133:          U_t = F(t,U) 
1134: .ve

1136:    Level: beginner

1138: .keywords: TS, problem type
1139: .seealso: TSSetUp(), TSProblemType, TS
1140: @*/
1141: PetscErrorCode  TSSetProblemType(TS ts, TSProblemType type)
1142: {

1147:   ts->problem_type = type;
1148:   if (type == TS_LINEAR) {
1149:     SNES snes;
1150:     TSGetSNES(ts,&snes);
1151:     SNESSetType(snes,SNESKSPONLY);
1152:   }
1153:   return(0);
1154: }

1158: /*@C
1159:   TSGetProblemType - Gets the type of problem to be solved.

1161:   Not collective

1163:   Input Parameter:
1164: . ts   - The TS

1166:   Output Parameter:
1167: . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1168: .vb
1169:          M U_t = A U
1170:          M(t) U_t = A(t) U
1171:          U_t = F(t,U)
1172: .ve

1174:    Level: beginner

1176: .keywords: TS, problem type
1177: .seealso: TSSetUp(), TSProblemType, TS
1178: @*/
1179: PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
1180: {
1184:   *type = ts->problem_type;
1185:   return(0);
1186: }

1190: /*@
1191:    TSSetUp - Sets up the internal data structures for the later use
1192:    of a timestepper.  

1194:    Collective on TS

1196:    Input Parameter:
1197: .  ts - the TS context obtained from TSCreate()

1199:    Notes:
1200:    For basic use of the TS solvers the user need not explicitly call
1201:    TSSetUp(), since these actions will automatically occur during
1202:    the call to TSStep().  However, if one wishes to control this
1203:    phase separately, TSSetUp() should be called after TSCreate()
1204:    and optional routines of the form TSSetXXX(), but before TSStep().  

1206:    Level: advanced

1208: .keywords: TS, timestep, setup

1210: .seealso: TSCreate(), TSStep(), TSDestroy()
1211: @*/
1212: PetscErrorCode  TSSetUp(TS ts)
1213: {

1218:   if (ts->setupcalled) return(0);

1220:   if (!((PetscObject)ts)->type_name) {
1221:     TSSetType(ts,TSEULER);
1222:   }
1223:   if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE;

1225:   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");

1227:   TSGetAdapt(ts,&ts->adapt);

1229:   if (ts->ops->setup) {
1230:     (*ts->ops->setup)(ts);
1231:   }

1233:   ts->setupcalled = PETSC_TRUE;
1234:   return(0);
1235: }

1239: /*@
1240:    TSReset - Resets a TS context and removes any allocated Vecs and Mats.

1242:    Collective on TS

1244:    Input Parameter:
1245: .  ts - the TS context obtained from TSCreate()

1247:    Level: beginner

1249: .keywords: TS, timestep, reset

1251: .seealso: TSCreate(), TSSetup(), TSDestroy()
1252: @*/
1253: PetscErrorCode  TSReset(TS ts)
1254: {

1259:   if (ts->ops->reset) {
1260:     (*ts->ops->reset)(ts);
1261:   }
1262:   if (ts->snes) {SNESReset(ts->snes);}
1263:   MatDestroy(&ts->Arhs);
1264:   MatDestroy(&ts->Brhs);
1265:   VecDestroy(&ts->Frhs);
1266:   VecDestroy(&ts->vec_sol);
1267:   VecDestroy(&ts->vatol);
1268:   VecDestroy(&ts->vrtol);
1269:   VecDestroyVecs(ts->nwork,&ts->work);
1270:   ts->setupcalled = PETSC_FALSE;
1271:   return(0);
1272: }

1276: /*@
1277:    TSDestroy - Destroys the timestepper context that was created
1278:    with TSCreate().

1280:    Collective on TS

1282:    Input Parameter:
1283: .  ts - the TS context obtained from TSCreate()

1285:    Level: beginner

1287: .keywords: TS, timestepper, destroy

1289: .seealso: TSCreate(), TSSetUp(), TSSolve()
1290: @*/
1291: PetscErrorCode  TSDestroy(TS *ts)
1292: {

1296:   if (!*ts) return(0);
1298:   if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; return(0);}

1300:   TSReset((*ts));

1302:   /* if memory was published with AMS then destroy it */
1303:   PetscObjectDepublish((*ts));
1304:   if ((*ts)->ops->destroy) {(*(*ts)->ops->destroy)((*ts));}

1306:   TSAdaptDestroy(&(*ts)->adapt);
1307:   SNESDestroy(&(*ts)->snes);
1308:   DMDestroy(&(*ts)->dm);
1309:   TSMonitorCancel((*ts));

1311:   PetscFree((*ts)->userops);

1313:   PetscHeaderDestroy(ts);
1314:   return(0);
1315: }

1319: /*@
1320:    TSGetSNES - Returns the SNES (nonlinear solver) associated with 
1321:    a TS (timestepper) context. Valid only for nonlinear problems.

1323:    Not Collective, but SNES is parallel if TS is parallel

1325:    Input Parameter:
1326: .  ts - the TS context obtained from TSCreate()

1328:    Output Parameter:
1329: .  snes - the nonlinear solver context

1331:    Notes:
1332:    The user can then directly manipulate the SNES context to set various
1333:    options, etc.  Likewise, the user can then extract and manipulate the
1334:    KSP, KSP, and PC contexts as well.

1336:    TSGetSNES() does not work for integrators that do not use SNES; in
1337:    this case TSGetSNES() returns PETSC_NULL in snes.

1339:    Level: beginner

1341: .keywords: timestep, get, SNES
1342: @*/
1343: PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
1344: {

1350:   if (!ts->snes) {
1351:     SNESCreate(((PetscObject)ts)->comm,&ts->snes);
1352:     PetscLogObjectParent(ts,ts->snes);
1353:     PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);
1354:     if (ts->dm) {SNESSetDM(ts->snes,ts->dm);}
1355:     if (ts->problem_type == TS_LINEAR) {
1356:       SNESSetType(ts->snes,SNESKSPONLY);
1357:     }
1358:   }
1359:   *snes = ts->snes;
1360:   return(0);
1361: }

1365: /*@
1366:    TSGetKSP - Returns the KSP (linear solver) associated with
1367:    a TS (timestepper) context.

1369:    Not Collective, but KSP is parallel if TS is parallel

1371:    Input Parameter:
1372: .  ts - the TS context obtained from TSCreate()

1374:    Output Parameter:
1375: .  ksp - the nonlinear solver context

1377:    Notes:
1378:    The user can then directly manipulate the KSP context to set various
1379:    options, etc.  Likewise, the user can then extract and manipulate the
1380:    KSP and PC contexts as well.

1382:    TSGetKSP() does not work for integrators that do not use KSP;
1383:    in this case TSGetKSP() returns PETSC_NULL in ksp.

1385:    Level: beginner

1387: .keywords: timestep, get, KSP
1388: @*/
1389: PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
1390: {
1392:   SNES           snes;

1397:   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
1398:   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1399:   TSGetSNES(ts,&snes);
1400:   SNESGetKSP(snes,ksp);
1401:   return(0);
1402: }

1404: /* ----------- Routines to set solver parameters ---------- */

1408: /*@
1409:    TSGetDuration - Gets the maximum number of timesteps to use and
1410:    maximum time for iteration.

1412:    Not Collective

1414:    Input Parameters:
1415: +  ts       - the TS context obtained from TSCreate()
1416: .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1417: -  maxtime  - final time to iterate to, or PETSC_NULL

1419:    Level: intermediate

1421: .keywords: TS, timestep, get, maximum, iterations, time
1422: @*/
1423: PetscErrorCode  TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1424: {
1427:   if (maxsteps) {
1429:     *maxsteps = ts->max_steps;
1430:   }
1431:   if (maxtime) {
1433:     *maxtime  = ts->max_time;
1434:   }
1435:   return(0);
1436: }

1440: /*@
1441:    TSSetDuration - Sets the maximum number of timesteps to use and
1442:    maximum time for iteration.

1444:    Logically Collective on TS

1446:    Input Parameters:
1447: +  ts - the TS context obtained from TSCreate()
1448: .  maxsteps - maximum number of iterations to use
1449: -  maxtime - final time to iterate to

1451:    Options Database Keys:
1452: .  -ts_max_steps <maxsteps> - Sets maxsteps
1453: .  -ts_final_time <maxtime> - Sets maxtime

1455:    Notes:
1456:    The default maximum number of iterations is 5000. Default time is 5.0

1458:    Level: intermediate

1460: .keywords: TS, timestep, set, maximum, iterations

1462: .seealso: TSSetExactFinalTime()
1463: @*/
1464: PetscErrorCode  TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1465: {
1470:   if (maxsteps >= 0) ts->max_steps = maxsteps;
1471:   if (maxtime != PETSC_DEFAULT) ts->max_time  = maxtime;
1472:   return(0);
1473: }

1477: /*@
1478:    TSSetSolution - Sets the initial solution vector
1479:    for use by the TS routines.

1481:    Logically Collective on TS and Vec

1483:    Input Parameters:
1484: +  ts - the TS context obtained from TSCreate()
1485: -  x - the solution vector

1487:    Level: beginner

1489: .keywords: TS, timestep, set, solution, initial conditions
1490: @*/
1491: PetscErrorCode  TSSetSolution(TS ts,Vec x)
1492: {
1494:   DM             dm;

1499:   PetscObjectReference((PetscObject)x);
1500:   VecDestroy(&ts->vec_sol);
1501:   ts->vec_sol = x;
1502:   TSGetDM(ts,&dm);
1503:   DMShellSetGlobalVector(dm,x);
1504:   return(0);
1505: }

1509: /*@C
1510:   TSSetPreStep - Sets the general-purpose function
1511:   called once at the beginning of each time step.

1513:   Logically Collective on TS

1515:   Input Parameters:
1516: + ts   - The TS context obtained from TSCreate()
1517: - func - The function

1519:   Calling sequence of func:
1520: . func (TS ts);

1522:   Level: intermediate

1524:   Note:
1525:   If a step is rejected, TSStep() will call this routine again before each attempt.
1526:   The last completed time step number can be queried using TSGetTimeStepNumber(), the
1527:   size of the step being attempted can be obtained using TSGetTimeStep().

1529: .keywords: TS, timestep
1530: .seealso: TSSetPreStage(), TSSetPostStep(), TSStep()
1531: @*/
1532: PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1533: {
1536:   ts->ops->prestep = func;
1537:   return(0);
1538: }

1542: /*@
1543:   TSPreStep - Runs the user-defined pre-step function.

1545:   Collective on TS

1547:   Input Parameters:
1548: . ts   - The TS context obtained from TSCreate()

1550:   Notes:
1551:   TSPreStep() is typically used within time stepping implementations,
1552:   so most users would not generally call this routine themselves.

1554:   Level: developer

1556: .keywords: TS, timestep
1557: .seealso: TSSetPreStep(), TSPreStage(), TSPostStep()
1558: @*/
1559: PetscErrorCode  TSPreStep(TS ts)
1560: {

1565:   if (ts->ops->prestep) {
1566:     PetscStackPush("TS PreStep function");
1567:     (*ts->ops->prestep)(ts);
1568:     PetscStackPop;
1569:   }
1570:   return(0);
1571: }

1575: /*@C
1576:   TSSetPreStage - Sets the general-purpose function
1577:   called once at the beginning of each stage.

1579:   Logically Collective on TS

1581:   Input Parameters:
1582: + ts   - The TS context obtained from TSCreate()
1583: - func - The function

1585:   Calling sequence of func:
1586: . PetscErrorCode func(TS ts, PetscReal stagetime);

1588:   Level: intermediate

1590:   Note:
1591:   There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
1592:   The time step number being computed can be queried using TSGetTimeStepNumber() and the total size of the step being
1593:   attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime().

1595: .keywords: TS, timestep
1596: .seealso: TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
1597: @*/
1598: PetscErrorCode  TSSetPreStage(TS ts, PetscErrorCode (*func)(TS,PetscReal))
1599: {
1602:   ts->ops->prestage = func;
1603:   return(0);
1604: }

1608: /*@
1609:   TSPreStage - Runs the user-defined pre-stage function set using TSSetPreStage()

1611:   Collective on TS

1613:   Input Parameters:
1614: . ts   - The TS context obtained from TSCreate()

1616:   Notes:
1617:   TSPreStage() is typically used within time stepping implementations,
1618:   most users would not generally call this routine themselves.

1620:   Level: developer

1622: .keywords: TS, timestep
1623: .seealso: TSSetPreStep(), TSPreStep(), TSPostStep()
1624: @*/
1625: PetscErrorCode  TSPreStage(TS ts, PetscReal stagetime)
1626: {

1631:   if (ts->ops->prestage) {
1632:     PetscStackPush("TS PreStage function");
1633:     (*ts->ops->prestage)(ts,stagetime);
1634:     PetscStackPop;
1635:   }
1636:   return(0);
1637: }

1641: /*@C
1642:   TSSetPostStep - Sets the general-purpose function
1643:   called once at the end of each time step.

1645:   Logically Collective on TS

1647:   Input Parameters:
1648: + ts   - The TS context obtained from TSCreate()
1649: - func - The function

1651:   Calling sequence of func:
1652: $ func (TS ts);

1654:   Level: intermediate

1656: .keywords: TS, timestep
1657: .seealso: TSSetPreStep(), TSSetPreStage(), TSGetTimeStep(), TSGetTimeStepNumber(), TSGetTime()
1658: @*/
1659: PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1660: {
1663:   ts->ops->poststep = func;
1664:   return(0);
1665: }

1669: /*@
1670:   TSPostStep - Runs the user-defined post-step function.

1672:   Collective on TS

1674:   Input Parameters:
1675: . ts   - The TS context obtained from TSCreate()

1677:   Notes:
1678:   TSPostStep() is typically used within time stepping implementations,
1679:   so most users would not generally call this routine themselves.

1681:   Level: developer

1683: .keywords: TS, timestep
1684: @*/
1685: PetscErrorCode  TSPostStep(TS ts)
1686: {

1691:   if (ts->ops->poststep) {
1692:     PetscStackPush("TS PostStep function");
1693:     (*ts->ops->poststep)(ts);
1694:     PetscStackPop;
1695:   }
1696:   return(0);
1697: }

1699: /* ------------ Routines to set performance monitoring options ----------- */

1703: /*@C
1704:    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1705:    timestep to display the iteration's  progress.

1707:    Logically Collective on TS

1709:    Input Parameters:
1710: +  ts - the TS context obtained from TSCreate()
1711: .  monitor - monitoring routine
1712: .  mctx - [optional] user-defined context for private data for the 
1713:              monitor routine (use PETSC_NULL if no context is desired)
1714: -  monitordestroy - [optional] routine that frees monitor context
1715:           (may be PETSC_NULL)

1717:    Calling sequence of monitor:
1718: $    int monitor(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)

1720: +    ts - the TS context
1721: .    steps - iteration number
1722: .    time - current time
1723: .    x - current iterate
1724: -    mctx - [optional] monitoring context

1726:    Notes:
1727:    This routine adds an additional monitor to the list of monitors that 
1728:    already has been loaded.

1730:    Fortran notes: Only a single monitor function can be set for each TS object

1732:    Level: intermediate

1734: .keywords: TS, timestep, set, monitor

1736: .seealso: TSMonitorDefault(), TSMonitorCancel()
1737: @*/
1738: PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
1739: {
1742:   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1743:   ts->monitor[ts->numbermonitors]           = monitor;
1744:   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1745:   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1746:   return(0);
1747: }

1751: /*@C
1752:    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.

1754:    Logically Collective on TS

1756:    Input Parameters:
1757: .  ts - the TS context obtained from TSCreate()

1759:    Notes:
1760:    There is no way to remove a single, specific monitor.

1762:    Level: intermediate

1764: .keywords: TS, timestep, set, monitor

1766: .seealso: TSMonitorDefault(), TSMonitorSet()
1767: @*/
1768: PetscErrorCode  TSMonitorCancel(TS ts)
1769: {
1771:   PetscInt       i;

1775:   for (i=0; i<ts->numbermonitors; i++) {
1776:     if (ts->mdestroy[i]) {
1777:       (*ts->mdestroy[i])(&ts->monitorcontext[i]);
1778:     }
1779:   }
1780:   ts->numbermonitors = 0;
1781:   return(0);
1782: }

1786: /*@
1787:    TSMonitorDefault - Sets the Default monitor

1789:    Level: intermediate

1791: .keywords: TS, set, monitor

1793: .seealso: TSMonitorDefault(), TSMonitorSet()
1794: @*/
1795: PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
1796: {
1798:   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);

1801:   PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);
1802:   PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g\n",step,(double)ts->time_step,(double)ptime);
1803:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);
1804:   return(0);
1805: }

1809: /*@
1810:    TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available.

1812:    Logically Collective on TS

1814:    Input Argument:
1815: .  ts - time stepping context

1817:    Output Argument:
1818: .  flg - PETSC_TRUE or PETSC_FALSE

1820:    Level: intermediate

1822: .keywords: TS, set

1824: .seealso: TSInterpolate(), TSSetPostStep()
1825: @*/
1826: PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg)
1827: {

1831:   ts->retain_stages = flg;
1832:   return(0);
1833: }

1837: /*@
1838:    TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval

1840:    Collective on TS

1842:    Input Argument:
1843: +  ts - time stepping context
1844: -  t - time to interpolate to

1846:    Output Argument:
1847: .  X - state at given time

1849:    Notes:
1850:    The user should call TSSetRetainStages() before taking a step in which interpolation will be requested.

1852:    Level: intermediate

1854:    Developer Notes:
1855:    TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints.

1857: .keywords: TS, set

1859: .seealso: TSSetRetainStages(), TSSetPostStep()
1860: @*/
1861: PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec X)
1862: {

1867:   if (t < ts->ptime - ts->time_step_prev || t > ts->ptime) SETERRQ3(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested time %G not in last time steps [%G,%G]",t,ts->ptime-ts->time_step_prev,ts->ptime);
1868:   if (!ts->ops->interpolate) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
1869:   (*ts->ops->interpolate)(ts,t,X);
1870:   return(0);
1871: }

1875: /*@
1876:    TSStep - Steps one time step

1878:    Collective on TS

1880:    Input Parameter:
1881: .  ts - the TS context obtained from TSCreate()

1883:    Level: intermediate

1885:    Notes:
1886:    The hook set using TSSetPreStep() is called before each attempt to take the step. In general, the time step size may
1887:    be changed due to adaptive error controller or solve failures. Note that steps may contain multiple stages.

1889: .keywords: TS, timestep, solve

1891: .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage()
1892: @*/
1893: PetscErrorCode  TSStep(TS ts)
1894: {
1895:   PetscReal      ptime_prev;

1900:   TSSetUp(ts);

1902:   ts->reason = TS_CONVERGED_ITERATING;

1904:   ptime_prev = ts->ptime;
1905:   PetscLogEventBegin(TS_Step,ts,0,0,0);
1906:   (*ts->ops->step)(ts);
1907:   PetscLogEventEnd(TS_Step,ts,0,0,0);
1908:   ts->time_step_prev = ts->ptime - ptime_prev;

1910:   if (ts->reason < 0) {
1911:     if (ts->errorifstepfailed) {
1912:       if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) {
1913:         SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_snes_failures or make negative to attempt recovery",TSConvergedReasons[ts->reason]);
1914:       } else SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
1915:     }
1916:   } else if (!ts->reason) {
1917:     if (ts->steps >= ts->max_steps)
1918:       ts->reason = TS_CONVERGED_ITS;
1919:     else if (ts->ptime >= ts->max_time)
1920:       ts->reason = TS_CONVERGED_TIME;
1921:   }

1923:   return(0);
1924: }

1928: /*@
1929:    TSEvaluateStep - Evaluate the solution at the end of a time step with a given order of accuracy.

1931:    Collective on TS

1933:    Input Arguments:
1934: +  ts - time stepping context
1935: .  order - desired order of accuracy
1936: -  done - whether the step was evaluated at this order (pass PETSC_NULL to generate an error if not available)

1938:    Output Arguments:
1939: .  X - state at the end of the current step

1941:    Level: advanced

1943:    Notes:
1944:    This function cannot be called until all stages have been evaluated.
1945:    It is normally called by adaptive controllers before a step has been accepted and may also be called by the user after TSStep() has returned.

1947: .seealso: TSStep(), TSAdapt
1948: @*/
1949: PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec X,PetscBool *done)
1950: {

1957:   if (!ts->ops->evaluatestep) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSEvaluateStep not implemented for type '%s'",((PetscObject)ts)->type_name);
1958:   (*ts->ops->evaluatestep)(ts,order,X,done);
1959:   return(0);
1960: }

1964: /*@
1965:    TSSolve - Steps the requested number of timesteps.

1967:    Collective on TS

1969:    Input Parameter:
1970: +  ts - the TS context obtained from TSCreate()
1971: -  x - the solution vector

1973:    Output Parameter:
1974: .  ftime - time of the state vector x upon completion

1976:    Level: beginner

1978:    Notes:
1979:    The final time returned by this function may be different from the time of the internally
1980:    held state accessible by TSGetSolution() and TSGetTime() because the method may have
1981:    stepped over the final time.

1983: .keywords: TS, timestep, solve

1985: .seealso: TSCreate(), TSSetSolution(), TSStep()
1986: @*/
1987: PetscErrorCode TSSolve(TS ts,Vec x,PetscReal *ftime)
1988: {
1989:   PetscBool      flg;
1990:   char           filename[PETSC_MAX_PATH_LEN];
1991:   PetscViewer    viewer;

1997:   if (ts->exact_final_time) {   /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */
1998:     if (!ts->vec_sol || x == ts->vec_sol) {
1999:       Vec y;
2000:       VecDuplicate(x,&y);
2001:       TSSetSolution(ts,y);
2002:       VecDestroy(&y); /* grant ownership */
2003:     }
2004:     VecCopy(x,ts->vec_sol);
2005:   } else {
2006:     TSSetSolution(ts,x);
2007:   }
2008:   TSSetUp(ts);
2009:   /* reset time step and iteration counters */
2010:   ts->steps = 0;
2011:   ts->ksp_its = 0;
2012:   ts->snes_its = 0;
2013:   ts->num_snes_failures = 0;
2014:   ts->reject = 0;
2015:   ts->reason = TS_CONVERGED_ITERATING;

2017:   if (ts->ops->solve) {         /* This private interface is transitional and should be removed when all implementations are updated. */
2018:     (*ts->ops->solve)(ts);
2019:     VecCopy(ts->vec_sol,x);
2020:     if (ftime) *ftime = ts->ptime;
2021:   } else {
2022:     /* steps the requested number of timesteps. */
2023:     TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
2024:     if (ts->steps >= ts->max_steps)
2025:       ts->reason = TS_CONVERGED_ITS;
2026:     else if (ts->ptime >= ts->max_time)
2027:       ts->reason = TS_CONVERGED_TIME;
2028:     while (!ts->reason) {
2029:       TSStep(ts);
2030:       TSPostStep(ts);
2031:       TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
2032:     }
2033:     if (ts->exact_final_time && ts->ptime > ts->max_time) {
2034:       TSInterpolate(ts,ts->max_time,x);
2035:       if (ftime) *ftime = ts->max_time;
2036:     } else {
2037:       VecCopy(ts->vec_sol,x);
2038:       if (ftime) *ftime = ts->ptime;
2039:     }
2040:   }
2041:   PetscOptionsGetString(((PetscObject)ts)->prefix,"-ts_view",filename,PETSC_MAX_PATH_LEN,&flg);
2042:   if (flg && !PetscPreLoadingOn) {
2043:     PetscViewerASCIIOpen(((PetscObject)ts)->comm,filename,&viewer);
2044:     TSView(ts,viewer);
2045:     PetscViewerDestroy(&viewer);
2046:   }
2047:   return(0);
2048: }

2052: /*@
2053:    TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet()

2055:    Collective on TS

2057:    Input Parameters:
2058: +  ts - time stepping context obtained from TSCreate()
2059: .  step - step number that has just completed
2060: .  ptime - model time of the state
2061: -  x - state at the current model time

2063:    Notes:
2064:    TSMonitor() is typically used within the time stepping implementations.
2065:    Users might call this function when using the TSStep() interface instead of TSSolve().

2067:    Level: advanced

2069: .keywords: TS, timestep
2070: @*/
2071: PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
2072: {
2074:   PetscInt       i,n = ts->numbermonitors;

2077:   for (i=0; i<n; i++) {
2078:     (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);
2079:   }
2080:   return(0);
2081: }

2083: /* ------------------------------------------------------------------------*/

2087: /*@C
2088:    TSMonitorLGCreate - Creates a line graph context for use with 
2089:    TS to monitor convergence of preconditioned residual norms.

2091:    Collective on TS

2093:    Input Parameters:
2094: +  host - the X display to open, or null for the local machine
2095: .  label - the title to put in the title bar
2096: .  x, y - the screen coordinates of the upper left coordinate of the window
2097: -  m, n - the screen width and height in pixels

2099:    Output Parameter:
2100: .  draw - the drawing context

2102:    Options Database Key:
2103: .  -ts_monitor_draw - automatically sets line graph monitor

2105:    Notes: 
2106:    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().

2108:    Level: intermediate

2110: .keywords: TS, monitor, line graph, residual, seealso

2112: .seealso: TSMonitorLGDestroy(), TSMonitorSet()

2114: @*/
2115: PetscErrorCode  TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
2116: {
2117:   PetscDraw      win;

2121:   PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);
2122:   PetscDrawSetType(win,PETSC_DRAW_X);
2123:   PetscDrawLGCreate(win,1,draw);
2124:   PetscDrawLGIndicateDataPoints(*draw);

2126:   PetscLogObjectParent(*draw,win);
2127:   return(0);
2128: }

2132: PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
2133: {
2134:   PetscDrawLG    lg = (PetscDrawLG) monctx;
2135:   PetscReal      x,y = ptime;

2139:   if (!monctx) {
2140:     MPI_Comm    comm;
2141:     PetscViewer viewer;

2143:     PetscObjectGetComm((PetscObject)ts,&comm);
2144:     viewer = PETSC_VIEWER_DRAW_(comm);
2145:     PetscViewerDrawGetDrawLG(viewer,0,&lg);
2146:   }

2148:   if (!n) {PetscDrawLGReset(lg);}
2149:   x = (PetscReal)n;
2150:   PetscDrawLGAddPoint(lg,&x,&y);
2151:   if (n < 20 || (n % 5)) {
2152:     PetscDrawLGDraw(lg);
2153:   }
2154:   return(0);
2155: }

2159: /*@C
2160:    TSMonitorLGDestroy - Destroys a line graph context that was created 
2161:    with TSMonitorLGCreate().

2163:    Collective on PetscDrawLG

2165:    Input Parameter:
2166: .  draw - the drawing context

2168:    Level: intermediate

2170: .keywords: TS, monitor, line graph, destroy

2172: .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
2173: @*/
2174: PetscErrorCode  TSMonitorLGDestroy(PetscDrawLG *drawlg)
2175: {
2176:   PetscDraw      draw;

2180:   PetscDrawLGGetDraw(*drawlg,&draw);
2181:   PetscDrawDestroy(&draw);
2182:   PetscDrawLGDestroy(drawlg);
2183:   return(0);
2184: }

2188: /*@
2189:    TSGetTime - Gets the time of the most recently completed step.

2191:    Not Collective

2193:    Input Parameter:
2194: .  ts - the TS context obtained from TSCreate()

2196:    Output Parameter:
2197: .  t  - the current time

2199:    Level: beginner

2201:    Note:
2202:    When called during time step evaluation (e.g. during residual evaluation or via hooks set using TSSetPreStep(),
2203:    TSSetPreStage(), or TSSetPostStep()), the time is the time at the start of the step being evaluated.

2205: .seealso: TSSetInitialTimeStep(), TSGetTimeStep()

2207: .keywords: TS, get, time
2208: @*/
2209: PetscErrorCode  TSGetTime(TS ts,PetscReal* t)
2210: {
2214:   *t = ts->ptime;
2215:   return(0);
2216: }

2220: /*@
2221:    TSSetTime - Allows one to reset the time.

2223:    Logically Collective on TS

2225:    Input Parameters:
2226: +  ts - the TS context obtained from TSCreate()
2227: -  time - the time

2229:    Level: intermediate

2231: .seealso: TSGetTime(), TSSetDuration()

2233: .keywords: TS, set, time
2234: @*/
2235: PetscErrorCode  TSSetTime(TS ts, PetscReal t)
2236: {
2240:   ts->ptime = t;
2241:   return(0);
2242: }

2246: /*@C
2247:    TSSetOptionsPrefix - Sets the prefix used for searching for all
2248:    TS options in the database.

2250:    Logically Collective on TS

2252:    Input Parameter:
2253: +  ts     - The TS context
2254: -  prefix - The prefix to prepend to all option names

2256:    Notes:
2257:    A hyphen (-) must NOT be given at the beginning of the prefix name.
2258:    The first character of all runtime options is AUTOMATICALLY the
2259:    hyphen.

2261:    Level: advanced

2263: .keywords: TS, set, options, prefix, database

2265: .seealso: TSSetFromOptions()

2267: @*/
2268: PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
2269: {
2271:   SNES           snes;

2275:   PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);
2276:   TSGetSNES(ts,&snes);
2277:   SNESSetOptionsPrefix(snes,prefix);
2278:   return(0);
2279: }


2284: /*@C
2285:    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
2286:    TS options in the database.

2288:    Logically Collective on TS

2290:    Input Parameter:
2291: +  ts     - The TS context
2292: -  prefix - The prefix to prepend to all option names

2294:    Notes:
2295:    A hyphen (-) must NOT be given at the beginning of the prefix name.
2296:    The first character of all runtime options is AUTOMATICALLY the
2297:    hyphen.

2299:    Level: advanced

2301: .keywords: TS, append, options, prefix, database

2303: .seealso: TSGetOptionsPrefix()

2305: @*/
2306: PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
2307: {
2309:   SNES           snes;

2313:   PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);
2314:   TSGetSNES(ts,&snes);
2315:   SNESAppendOptionsPrefix(snes,prefix);
2316:   return(0);
2317: }

2321: /*@C
2322:    TSGetOptionsPrefix - Sets the prefix used for searching for all
2323:    TS options in the database.

2325:    Not Collective

2327:    Input Parameter:
2328: .  ts - The TS context

2330:    Output Parameter:
2331: .  prefix - A pointer to the prefix string used

2333:    Notes: On the fortran side, the user should pass in a string 'prifix' of
2334:    sufficient length to hold the prefix.

2336:    Level: intermediate

2338: .keywords: TS, get, options, prefix, database

2340: .seealso: TSAppendOptionsPrefix()
2341: @*/
2342: PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
2343: {

2349:   PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);
2350:   return(0);
2351: }

2355: /*@C
2356:    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.

2358:    Not Collective, but parallel objects are returned if TS is parallel

2360:    Input Parameter:
2361: .  ts  - The TS context obtained from TSCreate()

2363:    Output Parameters:
2364: +  J   - The Jacobian J of F, where U_t = F(U,t)
2365: .  M   - The preconditioner matrix, usually the same as J
2366: .  func - Function to compute the Jacobian of the RHS
2367: -  ctx - User-defined context for Jacobian evaluation routine

2369:    Notes: You can pass in PETSC_NULL for any return argument you do not need.

2371:    Level: intermediate

2373: .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()

2375: .keywords: TS, timestep, get, matrix, Jacobian
2376: @*/
2377: PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx)
2378: {
2380:   SNES           snes;

2383:   TSGetSNES(ts,&snes);
2384:   SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);
2385:   if (func) *func = ts->userops->rhsjacobian;
2386:   if (ctx) *ctx = ts->jacP;
2387:   return(0);
2388: }

2392: /*@C
2393:    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.

2395:    Not Collective, but parallel objects are returned if TS is parallel

2397:    Input Parameter:
2398: .  ts  - The TS context obtained from TSCreate()

2400:    Output Parameters:
2401: +  A   - The Jacobian of F(t,U,U_t)
2402: .  B   - The preconditioner matrix, often the same as A
2403: .  f   - The function to compute the matrices
2404: - ctx - User-defined context for Jacobian evaluation routine

2406:    Notes: You can pass in PETSC_NULL for any return argument you do not need.

2408:    Level: advanced

2410: .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()

2412: .keywords: TS, timestep, get, matrix, Jacobian
2413: @*/
2414: PetscErrorCode  TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
2415: {
2417:   SNES           snes;

2420:   TSGetSNES(ts,&snes);
2421:   SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);
2422:   if (f) *f = ts->userops->ijacobian;
2423:   if (ctx) *ctx = ts->jacP;
2424:   return(0);
2425: }

2427: typedef struct {
2428:   PetscViewer viewer;
2429:   Vec         initialsolution;
2430:   PetscBool   showinitial;
2431: } TSMonitorSolutionCtx;

2435: /*@C
2436:    TSMonitorSolution - Monitors progress of the TS solvers by calling 
2437:    VecView() for the solution at each timestep

2439:    Collective on TS

2441:    Input Parameters:
2442: +  ts - the TS context
2443: .  step - current time-step
2444: .  ptime - current time
2445: -  dummy - either a viewer or PETSC_NULL

2447:    Level: intermediate

2449: .keywords: TS,  vector, monitor, view

2451: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
2452: @*/
2453: PetscErrorCode  TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
2454: {
2455:   PetscErrorCode       ierr;
2456:   TSMonitorSolutionCtx *ictx = (TSMonitorSolutionCtx*)dummy;

2459:   if (!step && ictx->showinitial) {
2460:     if (!ictx->initialsolution) {
2461:       VecDuplicate(x,&ictx->initialsolution);
2462:     }
2463:     VecCopy(x,ictx->initialsolution);
2464:   }
2465:   if (ictx->showinitial) {
2466:     PetscReal pause;
2467:     PetscViewerDrawGetPause(ictx->viewer,&pause);
2468:     PetscViewerDrawSetPause(ictx->viewer,0.0);
2469:     VecView(ictx->initialsolution,ictx->viewer);
2470:     PetscViewerDrawSetPause(ictx->viewer,pause);
2471:     PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);
2472:   }
2473:   VecView(x,ictx->viewer);
2474:   if (ictx->showinitial) {
2475:     PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);
2476:   }
2477:   return(0);
2478: }


2483: /*@C
2484:    TSMonitorSolutionDestroy - Destroys the monitor context for TSMonitorSolution

2486:    Collective on TS

2488:    Input Parameters:
2489: .    ctx - the monitor context

2491:    Level: intermediate

2493: .keywords: TS,  vector, monitor, view

2495: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
2496: @*/
2497: PetscErrorCode  TSMonitorSolutionDestroy(void **ctx)
2498: {
2499:   PetscErrorCode       ierr;
2500:   TSMonitorSolutionCtx *ictx = *(TSMonitorSolutionCtx**)ctx;
2501: 
2503:   PetscViewerDestroy(&ictx->viewer);
2504:   VecDestroy(&ictx->initialsolution);
2505:   PetscFree(ictx);
2506:   return(0);
2507: }

2511: /*@C
2512:    TSMonitorSolutionCreate - Creates the monitor context for TSMonitorSolution

2514:    Collective on TS

2516:    Input Parameter:
2517: .    ts - time-step context

2519:    Output Patameter:
2520: .    ctx - the monitor context

2522:    Level: intermediate

2524: .keywords: TS,  vector, monitor, view

2526: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
2527: @*/
2528: PetscErrorCode  TSMonitorSolutionCreate(TS ts,PetscViewer viewer,void **ctx)
2529: {
2530:   PetscErrorCode       ierr;
2531:   TSMonitorSolutionCtx *ictx;
2532: 
2534:   PetscNew(TSMonitorSolutionCtx,&ictx);
2535:   *ctx = (void*)ictx;
2536:   if (!viewer) {
2537:     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
2538:   }
2539:   PetscObjectReference((PetscObject)viewer);
2540:   ictx->viewer      = viewer;
2541:   ictx->showinitial = PETSC_FALSE;
2542:   PetscOptionsGetBool(((PetscObject)ts)->prefix,"-ts_monitor_solution_initial",&ictx->showinitial,PETSC_NULL);
2543:   return(0);
2544: }

2548: /*@
2549:    TSSetDM - Sets the DM that may be used by some preconditioners

2551:    Logically Collective on TS and DM

2553:    Input Parameters:
2554: +  ts - the preconditioner context
2555: -  dm - the dm

2557:    Level: intermediate


2560: .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
2561: @*/
2562: PetscErrorCode  TSSetDM(TS ts,DM dm)
2563: {
2565:   SNES           snes;

2569:   PetscObjectReference((PetscObject)dm);
2570:   DMDestroy(&ts->dm);
2571:   ts->dm = dm;
2572:   TSGetSNES(ts,&snes);
2573:   SNESSetDM(snes,dm);
2574:   return(0);
2575: }

2579: /*@
2580:    TSGetDM - Gets the DM that may be used by some preconditioners

2582:    Not Collective

2584:    Input Parameter:
2585: . ts - the preconditioner context

2587:    Output Parameter:
2588: .  dm - the dm

2590:    Level: intermediate


2593: .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
2594: @*/
2595: PetscErrorCode  TSGetDM(TS ts,DM *dm)
2596: {

2601:   if (!ts->dm) {
2602:     DMShellCreate(((PetscObject)ts)->comm,&ts->dm);
2603:     if (ts->snes) {SNESSetDM(ts->snes,ts->dm);}
2604:   }
2605:   *dm = ts->dm;
2606:   return(0);
2607: }

2611: /*@
2612:    SNESTSFormFunction - Function to evaluate nonlinear residual

2614:    Logically Collective on SNES

2616:    Input Parameter:
2617: + snes - nonlinear solver
2618: . X - the current state at which to evaluate the residual
2619: - ctx - user context, must be a TS

2621:    Output Parameter:
2622: . F - the nonlinear residual

2624:    Notes:
2625:    This function is not normally called by users and is automatically registered with the SNES used by TS.
2626:    It is most frequently passed to MatFDColoringSetFunction().

2628:    Level: advanced

2630: .seealso: SNESSetFunction(), MatFDColoringSetFunction()
2631: @*/
2632: PetscErrorCode  SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
2633: {
2634:   TS ts = (TS)ctx;

2642:   (ts->ops->snesfunction)(snes,X,F,ts);
2643:   return(0);
2644: }

2648: /*@
2649:    SNESTSFormJacobian - Function to evaluate the Jacobian

2651:    Collective on SNES

2653:    Input Parameter:
2654: + snes - nonlinear solver
2655: . X - the current state at which to evaluate the residual
2656: - ctx - user context, must be a TS

2658:    Output Parameter:
2659: + A - the Jacobian
2660: . B - the preconditioning matrix (may be the same as A)
2661: - flag - indicates any structure change in the matrix

2663:    Notes:
2664:    This function is not normally called by users and is automatically registered with the SNES used by TS.

2666:    Level: developer

2668: .seealso: SNESSetJacobian()
2669: @*/
2670: PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
2671: {
2672:   TS ts = (TS)ctx;

2684:   (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);
2685:   return(0);
2686: }

2690: /*@C
2691:    TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only

2693:    Collective on TS

2695:    Input Arguments:
2696: +  ts - time stepping context
2697: .  t - time at which to evaluate
2698: .  X - state at which to evaluate
2699: -  ctx - context

2701:    Output Arguments:
2702: .  F - right hand side

2704:    Level: intermediate

2706:    Notes:
2707:    This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems.
2708:    The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian().

2710: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
2711: @*/
2712: PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
2713: {
2715:   Mat Arhs,Brhs;
2716:   MatStructure flg2;

2719:   TSGetRHSMats_Private(ts,&Arhs,&Brhs);
2720:   TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);
2721:   MatMult(Arhs,X,F);
2722:   return(0);
2723: }

2727: /*@C
2728:    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.

2730:    Collective on TS

2732:    Input Arguments:
2733: +  ts - time stepping context
2734: .  t - time at which to evaluate
2735: .  X - state at which to evaluate
2736: -  ctx - context

2738:    Output Arguments:
2739: +  A - pointer to operator
2740: .  B - pointer to preconditioning matrix
2741: -  flg - matrix structure flag

2743:    Level: intermediate

2745:    Notes:
2746:    This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems.

2748: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
2749: @*/
2750: PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg,void *ctx)
2751: {

2754:   *flg = SAME_PRECONDITIONER;
2755:   return(0);
2756: }

2760: /*@C
2761:    TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only

2763:    Collective on TS

2765:    Input Arguments:
2766: +  ts - time stepping context
2767: .  t - time at which to evaluate
2768: .  X - state at which to evaluate
2769: .  Xdot - time derivative of state vector
2770: -  ctx - context

2772:    Output Arguments:
2773: .  F - left hand side

2775:    Level: intermediate

2777:    Notes:
2778:    The assumption here is that the left hand side is of the form A*Xdot (and not A*Xdot + B*X). For other cases, the
2779:    user is required to write their own TSComputeIFunction.
2780:    This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems.
2781:    The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian().

2783: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant()
2784: @*/
2785: PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
2786: {
2788:   Mat A,B;
2789:   MatStructure flg2;

2792:   TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);
2793:   TSComputeIJacobian(ts,t,X,Xdot,1.0,&A,&B,&flg2,PETSC_TRUE);
2794:   MatMult(A,Xdot,F);
2795:   return(0);
2796: }

2800: /*@C
2801:    TSComputeIJacobianConstant - Reuses a Jacobian that is time-independent.

2803:    Collective on TS

2805:    Input Arguments:
2806: +  ts - time stepping context
2807: .  t - time at which to evaluate
2808: .  X - state at which to evaluate
2809: .  Xdot - time derivative of state vector
2810: .  shift - shift to apply
2811: -  ctx - context

2813:    Output Arguments:
2814: +  A - pointer to operator
2815: .  B - pointer to preconditioning matrix
2816: -  flg - matrix structure flag

2818:    Level: intermediate

2820:    Notes:
2821:    This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems.

2823: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
2824: @*/
2825: PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx)
2826: {

2829:   *flg = SAME_PRECONDITIONER;
2830:   return(0);
2831: }


2836: /*@
2837:    TSGetConvergedReason - Gets the reason the TS iteration was stopped.

2839:    Not Collective

2841:    Input Parameter:
2842: .  ts - the TS context

2844:    Output Parameter:
2845: .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the 
2846:             manual pages for the individual convergence tests for complete lists

2848:    Level: intermediate

2850:    Notes:
2851:    Can only be called after the call to TSSolve() is complete.

2853: .keywords: TS, nonlinear, set, convergence, test

2855: .seealso: TSSetConvergenceTest(), TSConvergedReason
2856: @*/
2857: PetscErrorCode  TSGetConvergedReason(TS ts,TSConvergedReason *reason)
2858: {
2862:   *reason = ts->reason;
2863:   return(0);
2864: }

2868: /*@
2869:    TSGetSNESIterations - Gets the total number of nonlinear iterations
2870:    used by the time integrator.

2872:    Not Collective

2874:    Input Parameter:
2875: .  ts - TS context

2877:    Output Parameter:
2878: .  nits - number of nonlinear iterations

2880:    Notes:
2881:    This counter is reset to zero for each successive call to TSSolve().

2883:    Level: intermediate

2885: .keywords: TS, get, number, nonlinear, iterations

2887: .seealso:  TSGetKSPIterations()
2888: @*/
2889: PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits)
2890: {
2894:   *nits = ts->snes_its;
2895:   return(0);
2896: }

2900: /*@
2901:    TSGetKSPIterations - Gets the total number of linear iterations
2902:    used by the time integrator.

2904:    Not Collective

2906:    Input Parameter:
2907: .  ts - TS context

2909:    Output Parameter:
2910: .  lits - number of linear iterations

2912:    Notes:
2913:    This counter is reset to zero for each successive call to TSSolve().

2915:    Level: intermediate

2917: .keywords: TS, get, number, linear, iterations

2919: .seealso:  TSGetSNESIterations(), SNESGetKSPIterations()
2920: @*/
2921: PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits)
2922: {
2926:   *lits = ts->ksp_its;
2927:   return(0);
2928: }

2932: /*@
2933:    TSGetStepRejections - Gets the total number of rejected steps.

2935:    Not Collective

2937:    Input Parameter:
2938: .  ts - TS context

2940:    Output Parameter:
2941: .  rejects - number of steps rejected

2943:    Notes:
2944:    This counter is reset to zero for each successive call to TSSolve().

2946:    Level: intermediate

2948: .keywords: TS, get, number

2950: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails()
2951: @*/
2952: PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects)
2953: {
2957:   *rejects = ts->reject;
2958:   return(0);
2959: }

2963: /*@
2964:    TSGetSNESFailures - Gets the total number of failed SNES solves

2966:    Not Collective

2968:    Input Parameter:
2969: .  ts - TS context

2971:    Output Parameter:
2972: .  fails - number of failed nonlinear solves

2974:    Notes:
2975:    This counter is reset to zero for each successive call to TSSolve().

2977:    Level: intermediate

2979: .keywords: TS, get, number

2981: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures()
2982: @*/
2983: PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails)
2984: {
2988:   *fails = ts->num_snes_failures;
2989:   return(0);
2990: }

2994: /*@
2995:    TSSetMaxStepRejections - Sets the maximum number of step rejections before a step fails

2997:    Not Collective

2999:    Input Parameter:
3000: +  ts - TS context
3001: -  rejects - maximum number of rejected steps, pass -1 for unlimited

3003:    Notes:
3004:    The counter is reset to zero for each step

3006:    Options Database Key:
3007:  .  -ts_max_reject - Maximum number of step rejections before a step fails

3009:    Level: intermediate

3011: .keywords: TS, set, maximum, number

3013: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
3014: @*/
3015: PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects)
3016: {
3019:   ts->max_reject = rejects;
3020:   return(0);
3021: }

3025: /*@
3026:    TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves

3028:    Not Collective

3030:    Input Parameter:
3031: +  ts - TS context
3032: -  fails - maximum number of failed nonlinear solves, pass -1 for unlimited

3034:    Notes:
3035:    The counter is reset to zero for each successive call to TSSolve().

3037:    Options Database Key:
3038:  .  -ts_max_snes_failures - Maximum number of nonlinear solve failures

3040:    Level: intermediate

3042: .keywords: TS, set, maximum, number

3044: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason()
3045: @*/
3046: PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails)
3047: {
3050:   ts->max_snes_failures = fails;
3051:   return(0);
3052: }

3056: /*@
3057:    TSSetErrorIfStepFails - Error if no step succeeds

3059:    Not Collective

3061:    Input Parameter:
3062: +  ts - TS context
3063: -  err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure

3065:    Options Database Key:
3066:  .  -ts_error_if_step_fails - Error if no step succeeds

3068:    Level: intermediate

3070: .keywords: TS, set, error

3072: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
3073: @*/
3074: PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err)
3075: {
3078:   ts->errorifstepfailed = err;
3079:   return(0);
3080: }

3084: /*@C
3085:    TSMonitorSolutionBinary - Monitors progress of the TS solvers by VecView() for the solution at each timestep. Normally the viewer is a binary file

3087:    Collective on TS

3089:    Input Parameters:
3090: +  ts - the TS context
3091: .  step - current time-step
3092: .  ptime - current time
3093: .  x - current state
3094: -  viewer - binary viewer

3096:    Level: intermediate

3098: .keywords: TS,  vector, monitor, view

3100: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3101: @*/
3102: PetscErrorCode  TSMonitorSolutionBinary(TS ts,PetscInt step,PetscReal ptime,Vec x,void *viewer)
3103: {
3104:   PetscErrorCode       ierr;
3105:   PetscViewer          v = (PetscViewer)viewer;

3108:   VecView(x,v);
3109:   return(0);
3110: }

3114: /*@C
3115:    TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep.

3117:    Collective on TS

3119:    Input Parameters:
3120: +  ts - the TS context
3121: .  step - current time-step
3122: .  ptime - current time
3123: .  x - current state
3124: -  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)

3126:    Level: intermediate

3128:    Notes:
3129:    The VTK format does not allow writing multiple time steps in the same file, therefore a different file will be written for each time step.
3130:    These are named according to the file name template.

3132:    This function is normally passed as an argument to TSMonitorSet() along with TSMonitorSolutionVTKDestroy().

3134: .keywords: TS,  vector, monitor, view

3136: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3137: @*/
3138: PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec x,void *filenametemplate)
3139: {
3141:   char           filename[PETSC_MAX_PATH_LEN];
3142:   PetscViewer    viewer;

3145:   PetscSNPrintf(filename,sizeof filename,(const char*)filenametemplate,step);
3146:   PetscViewerVTKOpen(((PetscObject)ts)->comm,filename,FILE_MODE_WRITE,&viewer);
3147:   VecView(x,viewer);
3148:   PetscViewerDestroy(&viewer);
3149:   return(0);
3150: }

3154: /*@C
3155:    TSMonitorSolutionVTKDestroy - Destroy context for monitoring

3157:    Collective on TS

3159:    Input Parameters:
3160: .  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)

3162:    Level: intermediate

3164:    Note:
3165:    This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK().

3167: .keywords: TS,  vector, monitor, view

3169: .seealso: TSMonitorSet(), TSMonitorSolutionVTK()
3170: @*/
3171: PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
3172: {

3176:   PetscFree(*(char**)filenametemplate);
3177:   return(0);
3178: }

3182: /*@
3183:    TSGetAdapt - Get the adaptive controller context for the current method

3185:    Collective on TS if controller has not been created yet

3187:    Input Arguments:
3188: .  ts - time stepping context

3190:    Output Arguments:
3191: .  adapt - adaptive controller

3193:    Level: intermediate

3195: .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose()
3196: @*/
3197: PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt)
3198: {

3204:   if (!ts->adapt) {
3205:     TSAdaptCreate(((PetscObject)ts)->comm,&ts->adapt);
3206:     PetscLogObjectParent(ts,ts->adapt);
3207:     PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);
3208:   }
3209:   *adapt = ts->adapt;
3210:   return(0);
3211: }

3215: /*@
3216:    TSSetTolerances - Set tolerances for local truncation error when using adaptive controller

3218:    Logically Collective

3220:    Input Arguments:
3221: +  ts - time integration context
3222: .  atol - scalar absolute tolerances, PETSC_DECIDE to leave current value
3223: .  vatol - vector of absolute tolerances or PETSC_NULL, used in preference to atol if present
3224: .  rtol - scalar relative tolerances, PETSC_DECIDE to leave current value
3225: -  vrtol - vector of relative tolerances or PETSC_NULL, used in preference to atol if present

3227:    Level: beginner

3229: .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances()
3230: @*/
3231: PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol)
3232: {

3236:   if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol;
3237:   if (vatol) {
3238:     PetscObjectReference((PetscObject)vatol);
3239:     VecDestroy(&ts->vatol);
3240:     ts->vatol = vatol;
3241:   }
3242:   if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol;
3243:   if (vrtol) {
3244:     PetscObjectReference((PetscObject)vrtol);
3245:     VecDestroy(&ts->vrtol);
3246:     ts->vrtol = vrtol;
3247:   }
3248:   return(0);
3249: }

3253: /*@
3254:    TSGetTolerances - Get tolerances for local truncation error when using adaptive controller

3256:    Logically Collective

3258:    Input Arguments:
3259: .  ts - time integration context

3261:    Output Arguments:
3262: +  atol - scalar absolute tolerances, PETSC_NULL to ignore
3263: .  vatol - vector of absolute tolerances, PETSC_NULL to ignore
3264: .  rtol - scalar relative tolerances, PETSC_NULL to ignore
3265: -  vrtol - vector of relative tolerances, PETSC_NULL to ignore

3267:    Level: beginner

3269: .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances()
3270: @*/
3271: PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol)
3272: {

3275:   if (atol)  *atol  = ts->atol;
3276:   if (vatol) *vatol = ts->vatol;
3277:   if (rtol)  *rtol  = ts->rtol;
3278:   if (vrtol) *vrtol = ts->vrtol;
3279:   return(0);
3280: }

3284: /*@
3285:    TSErrorNormWRMS - compute a weighted norm of the difference between a vector and the current state

3287:    Collective on TS

3289:    Input Arguments:
3290: +  ts - time stepping context
3291: -  Y - state vector to be compared to ts->vec_sol

3293:    Output Arguments:
3294: .  norm - weighted norm, a value of 1.0 is considered small

3296:    Level: developer

3298: .seealso: TSSetTolerances()
3299: @*/
3300: PetscErrorCode TSErrorNormWRMS(TS ts,Vec Y,PetscReal *norm)
3301: {
3303:   PetscInt i,n,N;
3304:   const PetscScalar *x,*y;
3305:   Vec X;
3306:   PetscReal sum,gsum;

3312:   X = ts->vec_sol;
3314:   if (X == Y) SETERRQ(((PetscObject)X)->comm,PETSC_ERR_ARG_IDN,"Y cannot be the TS solution vector");

3316:   VecGetSize(X,&N);
3317:   VecGetLocalSize(X,&n);
3318:   VecGetArrayRead(X,&x);
3319:   VecGetArrayRead(Y,&y);
3320:   sum = 0.;
3321:   if (ts->vatol && ts->vrtol) {
3322:     const PetscScalar *atol,*rtol;
3323:     VecGetArrayRead(ts->vatol,&atol);
3324:     VecGetArrayRead(ts->vrtol,&rtol);
3325:     for (i=0; i<n; i++) {
3326:       PetscReal tol = PetscRealPart(atol[i]) + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3327:       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3328:     }
3329:     VecRestoreArrayRead(ts->vatol,&atol);
3330:     VecRestoreArrayRead(ts->vrtol,&rtol);
3331:   } else if (ts->vatol) {       /* vector atol, scalar rtol */
3332:     const PetscScalar *atol;
3333:     VecGetArrayRead(ts->vatol,&atol);
3334:     for (i=0; i<n; i++) {
3335:       PetscReal tol = PetscRealPart(atol[i]) + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3336:       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3337:     }
3338:     VecRestoreArrayRead(ts->vatol,&atol);
3339:   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
3340:     const PetscScalar *rtol;
3341:     VecGetArrayRead(ts->vrtol,&rtol);
3342:     for (i=0; i<n; i++) {
3343:       PetscReal tol = ts->atol + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3344:       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3345:     }
3346:     VecRestoreArrayRead(ts->vrtol,&rtol);
3347:   } else {                      /* scalar atol, scalar rtol */
3348:     for (i=0; i<n; i++) {
3349:       PetscReal tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3350:       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3351:     }
3352:   }
3353:   VecRestoreArrayRead(X,&x);
3354:   VecRestoreArrayRead(Y,&y);

3356:   MPI_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,((PetscObject)ts)->comm);
3357:   *norm = PetscSqrtReal(gsum / N);
3358:   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
3359:   return(0);
3360: }

3364: /*@
3365:    TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler

3367:    Logically Collective on TS

3369:    Input Arguments:
3370: +  ts - time stepping context
3371: -  cfltime - maximum stable time step if using forward Euler (value can be different on each process)

3373:    Note:
3374:    After calling this function, the global CFL time can be obtained by calling TSGetCFLTime()

3376:    Level: intermediate

3378: .seealso: TSGetCFLTime(), TSADAPTCFL
3379: @*/
3380: PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime)
3381: {

3385:   ts->cfltime_local = cfltime;
3386:   ts->cfltime = -1.;
3387:   return(0);
3388: }

3392: /*@
3393:    TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler

3395:    Collective on TS

3397:    Input Arguments:
3398: .  ts - time stepping context

3400:    Output Arguments:
3401: .  cfltime - maximum stable time step for forward Euler

3403:    Level: advanced

3405: .seealso: TSSetCFLTimeLocal()
3406: @*/
3407: PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime)
3408: {

3412:   if (ts->cfltime < 0) {
3413:     MPI_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,((PetscObject)ts)->comm);
3414:   }
3415:   *cfltime = ts->cfltime;
3416:   return(0);
3417: }

3421: /*@
3422:    TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu

3424:    Input Parameters:
3425: .  ts   - the TS context.
3426: .  xl   - lower bound.
3427: .  xu   - upper bound.

3429:    Notes:
3430:    If this routine is not called then the lower and upper bounds are set to 
3431:    SNES_VI_NINF and SNES_VI_INF respectively during SNESSetUp().

3433:    Level: advanced

3435: @*/
3436: PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
3437: {
3439:   SNES           snes;

3442:   TSGetSNES(ts,&snes);
3443:   SNESVISetVariableBounds(snes,xl,xu);
3444:   return(0);
3445: }

3447: #if defined(PETSC_HAVE_MATLAB_ENGINE)
3448: #include <mex.h>

3450: typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;

3454: /*
3455:    TSComputeFunction_Matlab - Calls the function that has been set with
3456:                          TSSetFunctionMatlab().

3458:    Collective on TS

3460:    Input Parameters:
3461: +  snes - the TS context
3462: -  x - input vector

3464:    Output Parameter:
3465: .  y - function vector, as set by TSSetFunction()

3467:    Notes:
3468:    TSComputeFunction() is typically used within nonlinear solvers
3469:    implementations, so most users would not generally call this routine
3470:    themselves.

3472:    Level: developer

3474: .keywords: TS, nonlinear, compute, function

3476: .seealso: TSSetFunction(), TSGetFunction()
3477: */
3478: PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx)
3479: {
3480:   PetscErrorCode   ierr;
3481:   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3482:   int              nlhs = 1,nrhs = 7;
3483:   mxArray          *plhs[1],*prhs[7];
3484:   long long int    lx = 0,lxdot = 0,ly = 0,ls = 0;


3494:   PetscMemcpy(&ls,&snes,sizeof(snes));
3495:   PetscMemcpy(&lx,&x,sizeof(x));
3496:   PetscMemcpy(&lxdot,&xdot,sizeof(xdot));
3497:   PetscMemcpy(&ly,&y,sizeof(x));
3498:   prhs[0] =  mxCreateDoubleScalar((double)ls);
3499:   prhs[1] =  mxCreateDoubleScalar(time);
3500:   prhs[2] =  mxCreateDoubleScalar((double)lx);
3501:   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
3502:   prhs[4] =  mxCreateDoubleScalar((double)ly);
3503:   prhs[5] =  mxCreateString(sctx->funcname);
3504:   prhs[6] =  sctx->ctx;
3505:    mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");
3506:    mxGetScalar(plhs[0]);
3507:   mxDestroyArray(prhs[0]);
3508:   mxDestroyArray(prhs[1]);
3509:   mxDestroyArray(prhs[2]);
3510:   mxDestroyArray(prhs[3]);
3511:   mxDestroyArray(prhs[4]);
3512:   mxDestroyArray(prhs[5]);
3513:   mxDestroyArray(plhs[0]);
3514:   return(0);
3515: }


3520: /*
3521:    TSSetFunctionMatlab - Sets the function evaluation routine and function
3522:    vector for use by the TS routines in solving ODEs
3523:    equations from MATLAB. Here the function is a string containing the name of a MATLAB function

3525:    Logically Collective on TS

3527:    Input Parameters:
3528: +  ts - the TS context
3529: -  func - function evaluation routine

3531:    Calling sequence of func:
3532: $    func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx);

3534:    Level: beginner

3536: .keywords: TS, nonlinear, set, function

3538: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3539: */
3540: PetscErrorCode  TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx)
3541: {
3542:   PetscErrorCode  ierr;
3543:   TSMatlabContext *sctx;

3546:   /* currently sctx is memory bleed */
3547:   PetscMalloc(sizeof(TSMatlabContext),&sctx);
3548:   PetscStrallocpy(func,&sctx->funcname);
3549:   /*
3550:      This should work, but it doesn't
3551:   sctx->ctx = ctx;
3552:   mexMakeArrayPersistent(sctx->ctx);
3553:   */
3554:   sctx->ctx = mxDuplicateArray(ctx);
3555:   TSSetIFunction(ts,PETSC_NULL,TSComputeFunction_Matlab,sctx);
3556:   return(0);
3557: }

3561: /*
3562:    TSComputeJacobian_Matlab - Calls the function that has been set with
3563:                          TSSetJacobianMatlab().

3565:    Collective on TS

3567:    Input Parameters:
3568: +  ts - the TS context
3569: .  x - input vector
3570: .  A, B - the matrices
3571: -  ctx - user context

3573:    Output Parameter:
3574: .  flag - structure of the matrix

3576:    Level: developer

3578: .keywords: TS, nonlinear, compute, function

3580: .seealso: TSSetFunction(), TSGetFunction()
3581: @*/
3582: PetscErrorCode  TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx)
3583: {
3584:   PetscErrorCode  ierr;
3585:   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3586:   int             nlhs = 2,nrhs = 9;
3587:   mxArray         *plhs[2],*prhs[9];
3588:   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;


3594:   /* call Matlab function in ctx with arguments x and y */

3596:   PetscMemcpy(&ls,&ts,sizeof(ts));
3597:   PetscMemcpy(&lx,&x,sizeof(x));
3598:   PetscMemcpy(&lxdot,&xdot,sizeof(x));
3599:   PetscMemcpy(&lA,A,sizeof(x));
3600:   PetscMemcpy(&lB,B,sizeof(x));
3601:   prhs[0] =  mxCreateDoubleScalar((double)ls);
3602:   prhs[1] =  mxCreateDoubleScalar((double)time);
3603:   prhs[2] =  mxCreateDoubleScalar((double)lx);
3604:   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
3605:   prhs[4] =  mxCreateDoubleScalar((double)shift);
3606:   prhs[5] =  mxCreateDoubleScalar((double)lA);
3607:   prhs[6] =  mxCreateDoubleScalar((double)lB);
3608:   prhs[7] =  mxCreateString(sctx->funcname);
3609:   prhs[8] =  sctx->ctx;
3610:    mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");
3611:    mxGetScalar(plhs[0]);
3612:   *flag   =  (MatStructure) mxGetScalar(plhs[1]);
3613:   mxDestroyArray(prhs[0]);
3614:   mxDestroyArray(prhs[1]);
3615:   mxDestroyArray(prhs[2]);
3616:   mxDestroyArray(prhs[3]);
3617:   mxDestroyArray(prhs[4]);
3618:   mxDestroyArray(prhs[5]);
3619:   mxDestroyArray(prhs[6]);
3620:   mxDestroyArray(prhs[7]);
3621:   mxDestroyArray(plhs[0]);
3622:   mxDestroyArray(plhs[1]);
3623:   return(0);
3624: }


3629: /*
3630:    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
3631:    vector for use by the TS routines in solving ODEs from MATLAB. Here the function is a string containing the name of a MATLAB function

3633:    Logically Collective on TS

3635:    Input Parameters:
3636: +  ts - the TS context
3637: .  A,B - Jacobian matrices
3638: .  func - function evaluation routine
3639: -  ctx - user context

3641:    Calling sequence of func:
3642: $    flag = func (TS ts,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx);


3645:    Level: developer

3647: .keywords: TS, nonlinear, set, function

3649: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3650: */
3651: PetscErrorCode  TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx)
3652: {
3653:   PetscErrorCode    ierr;
3654:   TSMatlabContext *sctx;

3657:   /* currently sctx is memory bleed */
3658:   PetscMalloc(sizeof(TSMatlabContext),&sctx);
3659:   PetscStrallocpy(func,&sctx->funcname);
3660:   /*
3661:      This should work, but it doesn't
3662:   sctx->ctx = ctx;
3663:   mexMakeArrayPersistent(sctx->ctx);
3664:   */
3665:   sctx->ctx = mxDuplicateArray(ctx);
3666:   TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);
3667:   return(0);
3668: }

3672: /*
3673:    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().

3675:    Collective on TS

3677: .seealso: TSSetFunction(), TSGetFunction()
3678: @*/
3679: PetscErrorCode  TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec x, void *ctx)
3680: {
3681:   PetscErrorCode  ierr;
3682:   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3683:   int             nlhs = 1,nrhs = 6;
3684:   mxArray         *plhs[1],*prhs[6];
3685:   long long int   lx = 0,ls = 0;


3691:   PetscMemcpy(&ls,&ts,sizeof(ts));
3692:   PetscMemcpy(&lx,&x,sizeof(x));
3693:   prhs[0] =  mxCreateDoubleScalar((double)ls);
3694:   prhs[1] =  mxCreateDoubleScalar((double)it);
3695:   prhs[2] =  mxCreateDoubleScalar((double)time);
3696:   prhs[3] =  mxCreateDoubleScalar((double)lx);
3697:   prhs[4] =  mxCreateString(sctx->funcname);
3698:   prhs[5] =  sctx->ctx;
3699:    mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");
3700:    mxGetScalar(plhs[0]);
3701:   mxDestroyArray(prhs[0]);
3702:   mxDestroyArray(prhs[1]);
3703:   mxDestroyArray(prhs[2]);
3704:   mxDestroyArray(prhs[3]);
3705:   mxDestroyArray(prhs[4]);
3706:   mxDestroyArray(plhs[0]);
3707:   return(0);
3708: }


3713: /*
3714:    TSMonitorSetMatlab - Sets the monitor function from Matlab

3716:    Level: developer

3718: .keywords: TS, nonlinear, set, function

3720: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3721: */
3722: PetscErrorCode  TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx)
3723: {
3724:   PetscErrorCode    ierr;
3725:   TSMatlabContext *sctx;

3728:   /* currently sctx is memory bleed */
3729:   PetscMalloc(sizeof(TSMatlabContext),&sctx);
3730:   PetscStrallocpy(func,&sctx->funcname);
3731:   /*
3732:      This should work, but it doesn't
3733:   sctx->ctx = ctx;
3734:   mexMakeArrayPersistent(sctx->ctx);
3735:   */
3736:   sctx->ctx = mxDuplicateArray(ctx);
3737:   TSMonitorSet(ts,TSMonitor_Matlab,sctx,PETSC_NULL);
3738:   return(0);
3739: }
3740: #endif