Actual source code: ts.c

petsc-3.3-p2 2012-07-13
  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 current number of timesteps.

954:    Not Collective

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

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

962:    Level: intermediate

964: .keywords: TS, timestep, get, iteration, number
965: @*/
966: PetscErrorCode  TSGetTimeStepNumber(TS ts,PetscInt* iter)
967: {
971:   *iter = ts->steps;
972:   return(0);
973: }

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

981:    Logically Collective on TS

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

988:    Level: intermediate

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

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

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

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

1011:    Logically Collective on TS

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

1017:    Level: intermediate

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

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

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

1039:   Logically Collective on TS

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

1045:    Level: beginner

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

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

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

1064:    Not Collective

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

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

1072:    Level: intermediate

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

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

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

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

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

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

1103:    Level: intermediate

1105: .seealso: TSGetTimeStep()

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

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

1124:   Not collective

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

1135:    Level: beginner

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

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

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

1160:   Not collective

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

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

1173:    Level: beginner

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

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

1193:    Collective on TS

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

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

1205:    Level: advanced

1207: .keywords: TS, timestep, setup

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

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

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

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

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

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

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

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

1241:    Collective on TS

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

1246:    Level: beginner

1248: .keywords: TS, timestep, reset

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

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

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

1279:    Collective on TS

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

1284:    Level: beginner

1286: .keywords: TS, timestepper, destroy

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

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

1299:   TSReset((*ts));

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

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

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

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

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

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

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

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

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

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

1338:    Level: beginner

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

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

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

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

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

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

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

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

1384:    Level: beginner

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

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

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

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

1411:    Not Collective

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

1418:    Level: intermediate

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

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

1443:    Logically Collective on TS

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

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

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

1457:    Level: intermediate

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

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

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

1480:    Logically Collective on TS and Vec

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

1486:    Level: beginner

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

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

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

1512:   Logically Collective on TS

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

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

1521:   Level: intermediate

1523: .keywords: TS, timestep
1524: @*/
1525: PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1526: {
1529:   ts->ops->prestep = func;
1530:   return(0);
1531: }

1535: /*@
1536:   TSPreStep - Runs the user-defined pre-step function.

1538:   Collective on TS

1540:   Input Parameters:
1541: . ts   - The TS context obtained from TSCreate()

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

1547:   Level: developer

1549: .keywords: TS, timestep
1550: @*/
1551: PetscErrorCode  TSPreStep(TS ts)
1552: {

1557:   if (ts->ops->prestep) {
1558:     PetscStackPush("TS PreStep function");
1559:     (*ts->ops->prestep)(ts);
1560:     PetscStackPop;
1561:   }
1562:   return(0);
1563: }

1567: /*@C
1568:   TSSetPostStep - Sets the general-purpose function
1569:   called once at the end of each time step.

1571:   Logically Collective on TS

1573:   Input Parameters:
1574: + ts   - The TS context obtained from TSCreate()
1575: - func - The function

1577:   Calling sequence of func:
1578: . func (TS ts);

1580:   Level: intermediate

1582: .keywords: TS, timestep
1583: @*/
1584: PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1585: {
1588:   ts->ops->poststep = func;
1589:   return(0);
1590: }

1594: /*@
1595:   TSPostStep - Runs the user-defined post-step function.

1597:   Collective on TS

1599:   Input Parameters:
1600: . ts   - The TS context obtained from TSCreate()

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

1606:   Level: developer

1608: .keywords: TS, timestep
1609: @*/
1610: PetscErrorCode  TSPostStep(TS ts)
1611: {

1616:   if (ts->ops->poststep) {
1617:     PetscStackPush("TS PostStep function");
1618:     (*ts->ops->poststep)(ts);
1619:     PetscStackPop;
1620:   }
1621:   return(0);
1622: }

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

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

1632:    Logically Collective on TS

1634:    Input Parameters:
1635: +  ts - the TS context obtained from TSCreate()
1636: .  monitor - monitoring routine
1637: .  mctx - [optional] user-defined context for private data for the 
1638:              monitor routine (use PETSC_NULL if no context is desired)
1639: -  monitordestroy - [optional] routine that frees monitor context
1640:           (may be PETSC_NULL)

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

1645: +    ts - the TS context
1646: .    steps - iteration number
1647: .    time - current time
1648: .    x - current iterate
1649: -    mctx - [optional] monitoring context

1651:    Notes:
1652:    This routine adds an additional monitor to the list of monitors that 
1653:    already has been loaded.

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

1657:    Level: intermediate

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

1661: .seealso: TSMonitorDefault(), TSMonitorCancel()
1662: @*/
1663: PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
1664: {
1667:   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1668:   ts->monitor[ts->numbermonitors]           = monitor;
1669:   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1670:   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1671:   return(0);
1672: }

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

1679:    Logically Collective on TS

1681:    Input Parameters:
1682: .  ts - the TS context obtained from TSCreate()

1684:    Notes:
1685:    There is no way to remove a single, specific monitor.

1687:    Level: intermediate

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

1691: .seealso: TSMonitorDefault(), TSMonitorSet()
1692: @*/
1693: PetscErrorCode  TSMonitorCancel(TS ts)
1694: {
1696:   PetscInt       i;

1700:   for (i=0; i<ts->numbermonitors; i++) {
1701:     if (ts->mdestroy[i]) {
1702:       (*ts->mdestroy[i])(&ts->monitorcontext[i]);
1703:     }
1704:   }
1705:   ts->numbermonitors = 0;
1706:   return(0);
1707: }

1711: /*@
1712:    TSMonitorDefault - Sets the Default monitor

1714:    Level: intermediate

1716: .keywords: TS, set, monitor

1718: .seealso: TSMonitorDefault(), TSMonitorSet()
1719: @*/
1720: PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
1721: {
1723:   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);

1726:   PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);
1727:   PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g\n",step,(double)ts->time_step,(double)ptime);
1728:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);
1729:   return(0);
1730: }

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

1737:    Logically Collective on TS

1739:    Input Argument:
1740: .  ts - time stepping context

1742:    Output Argument:
1743: .  flg - PETSC_TRUE or PETSC_FALSE

1745:    Level: intermediate

1747: .keywords: TS, set

1749: .seealso: TSInterpolate(), TSSetPostStep()
1750: @*/
1751: PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg)
1752: {

1756:   ts->retain_stages = flg;
1757:   return(0);
1758: }

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

1765:    Collective on TS

1767:    Input Argument:
1768: +  ts - time stepping context
1769: -  t - time to interpolate to

1771:    Output Argument:
1772: .  X - state at given time

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

1777:    Level: intermediate

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

1782: .keywords: TS, set

1784: .seealso: TSSetRetainStages(), TSSetPostStep()
1785: @*/
1786: PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec X)
1787: {

1792:   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);
1793:   if (!ts->ops->interpolate) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
1794:   (*ts->ops->interpolate)(ts,t,X);
1795:   return(0);
1796: }

1800: /*@
1801:    TSStep - Steps one time step

1803:    Collective on TS

1805:    Input Parameter:
1806: .  ts - the TS context obtained from TSCreate()

1808:    Level: intermediate

1810: .keywords: TS, timestep, solve

1812: .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve()
1813: @*/
1814: PetscErrorCode  TSStep(TS ts)
1815: {
1816:   PetscReal      ptime_prev;

1821:   TSSetUp(ts);

1823:   ts->reason = TS_CONVERGED_ITERATING;

1825:   ptime_prev = ts->ptime;
1826:   PetscLogEventBegin(TS_Step,ts,0,0,0);
1827:   (*ts->ops->step)(ts);
1828:   PetscLogEventEnd(TS_Step,ts,0,0,0);
1829:   ts->time_step_prev = ts->ptime - ptime_prev;

1831:   if (ts->reason < 0) {
1832:     if (ts->errorifstepfailed) {
1833:       if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) {
1834:         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]);
1835:       } else SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
1836:     }
1837:   } else if (!ts->reason) {
1838:     if (ts->steps >= ts->max_steps)
1839:       ts->reason = TS_CONVERGED_ITS;
1840:     else if (ts->ptime >= ts->max_time)
1841:       ts->reason = TS_CONVERGED_TIME;
1842:   }

1844:   return(0);
1845: }

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

1852:    Collective on TS

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

1859:    Output Arguments:
1860: .  X - state at the end of the current step

1862:    Level: advanced

1864:    Notes:
1865:    This function cannot be called until all stages have been evaluated.
1866:    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.

1868: .seealso: TSStep(), TSAdapt
1869: @*/
1870: PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec X,PetscBool *done)
1871: {

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

1885: /*@
1886:    TSSolve - Steps the requested number of timesteps.

1888:    Collective on TS

1890:    Input Parameter:
1891: +  ts - the TS context obtained from TSCreate()
1892: -  x - the solution vector

1894:    Output Parameter:
1895: .  ftime - time of the state vector x upon completion

1897:    Level: beginner

1899:    Notes:
1900:    The final time returned by this function may be different from the time of the internally
1901:    held state accessible by TSGetSolution() and TSGetTime() because the method may have
1902:    stepped over the final time.

1904: .keywords: TS, timestep, solve

1906: .seealso: TSCreate(), TSSetSolution(), TSStep()
1907: @*/
1908: PetscErrorCode TSSolve(TS ts,Vec x,PetscReal *ftime)
1909: {
1910:   PetscBool      flg;
1911:   char           filename[PETSC_MAX_PATH_LEN];
1912:   PetscViewer    viewer;

1918:   if (ts->exact_final_time) {   /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */
1919:     if (!ts->vec_sol || x == ts->vec_sol) {
1920:       Vec y;
1921:       VecDuplicate(x,&y);
1922:       TSSetSolution(ts,y);
1923:       VecDestroy(&y); /* grant ownership */
1924:     }
1925:     VecCopy(x,ts->vec_sol);
1926:   } else {
1927:     TSSetSolution(ts,x);
1928:   }
1929:   TSSetUp(ts);
1930:   /* reset time step and iteration counters */
1931:   ts->steps = 0;
1932:   ts->ksp_its = 0;
1933:   ts->snes_its = 0;
1934:   ts->num_snes_failures = 0;
1935:   ts->reject = 0;
1936:   ts->reason = TS_CONVERGED_ITERATING;

1938:   if (ts->ops->solve) {         /* This private interface is transitional and should be removed when all implementations are updated. */
1939:     (*ts->ops->solve)(ts);
1940:     VecCopy(ts->vec_sol,x);
1941:     if (ftime) *ftime = ts->ptime;
1942:   } else {
1943:     /* steps the requested number of timesteps. */
1944:     TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
1945:     if (ts->steps >= ts->max_steps)
1946:       ts->reason = TS_CONVERGED_ITS;
1947:     else if (ts->ptime >= ts->max_time)
1948:       ts->reason = TS_CONVERGED_TIME;
1949:     while (!ts->reason) {
1950:       TSPreStep(ts);
1951:       TSStep(ts);
1952:       TSPostStep(ts);
1953:       TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
1954:     }
1955:     if (ts->exact_final_time && ts->ptime > ts->max_time) {
1956:       TSInterpolate(ts,ts->max_time,x);
1957:       if (ftime) *ftime = ts->max_time;
1958:     } else {
1959:       VecCopy(ts->vec_sol,x);
1960:       if (ftime) *ftime = ts->ptime;
1961:     }
1962:   }
1963:   PetscOptionsGetString(((PetscObject)ts)->prefix,"-ts_view",filename,PETSC_MAX_PATH_LEN,&flg);
1964:   if (flg && !PetscPreLoadingOn) {
1965:     PetscViewerASCIIOpen(((PetscObject)ts)->comm,filename,&viewer);
1966:     TSView(ts,viewer);
1967:     PetscViewerDestroy(&viewer);
1968:   }
1969:   return(0);
1970: }

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

1977:    Collective on TS

1979:    Input Parameters:
1980: +  ts - time stepping context obtained from TSCreate()
1981: .  step - step number that has just completed
1982: .  ptime - model time of the state
1983: -  x - state at the current model time

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

1989:    Level: advanced

1991: .keywords: TS, timestep
1992: @*/
1993: PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1994: {
1996:   PetscInt       i,n = ts->numbermonitors;

1999:   for (i=0; i<n; i++) {
2000:     (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);
2001:   }
2002:   return(0);
2003: }

2005: /* ------------------------------------------------------------------------*/

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

2013:    Collective on TS

2015:    Input Parameters:
2016: +  host - the X display to open, or null for the local machine
2017: .  label - the title to put in the title bar
2018: .  x, y - the screen coordinates of the upper left coordinate of the window
2019: -  m, n - the screen width and height in pixels

2021:    Output Parameter:
2022: .  draw - the drawing context

2024:    Options Database Key:
2025: .  -ts_monitor_draw - automatically sets line graph monitor

2027:    Notes: 
2028:    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().

2030:    Level: intermediate

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

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

2036: @*/
2037: PetscErrorCode  TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
2038: {
2039:   PetscDraw      win;

2043:   PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);
2044:   PetscDrawSetType(win,PETSC_DRAW_X);
2045:   PetscDrawLGCreate(win,1,draw);
2046:   PetscDrawLGIndicateDataPoints(*draw);

2048:   PetscLogObjectParent(*draw,win);
2049:   return(0);
2050: }

2054: PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
2055: {
2056:   PetscDrawLG    lg = (PetscDrawLG) monctx;
2057:   PetscReal      x,y = ptime;

2061:   if (!monctx) {
2062:     MPI_Comm    comm;
2063:     PetscViewer viewer;

2065:     PetscObjectGetComm((PetscObject)ts,&comm);
2066:     viewer = PETSC_VIEWER_DRAW_(comm);
2067:     PetscViewerDrawGetDrawLG(viewer,0,&lg);
2068:   }

2070:   if (!n) {PetscDrawLGReset(lg);}
2071:   x = (PetscReal)n;
2072:   PetscDrawLGAddPoint(lg,&x,&y);
2073:   if (n < 20 || (n % 5)) {
2074:     PetscDrawLGDraw(lg);
2075:   }
2076:   return(0);
2077: }

2081: /*@C
2082:    TSMonitorLGDestroy - Destroys a line graph context that was created 
2083:    with TSMonitorLGCreate().

2085:    Collective on PetscDrawLG

2087:    Input Parameter:
2088: .  draw - the drawing context

2090:    Level: intermediate

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

2094: .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
2095: @*/
2096: PetscErrorCode  TSMonitorLGDestroy(PetscDrawLG *drawlg)
2097: {
2098:   PetscDraw      draw;

2102:   PetscDrawLGGetDraw(*drawlg,&draw);
2103:   PetscDrawDestroy(&draw);
2104:   PetscDrawLGDestroy(drawlg);
2105:   return(0);
2106: }

2110: /*@
2111:    TSGetTime - Gets the current time.

2113:    Not Collective

2115:    Input Parameter:
2116: .  ts - the TS context obtained from TSCreate()

2118:    Output Parameter:
2119: .  t  - the current time

2121:    Level: beginner

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

2125: .keywords: TS, get, time
2126: @*/
2127: PetscErrorCode  TSGetTime(TS ts,PetscReal* t)
2128: {
2132:   *t = ts->ptime;
2133:   return(0);
2134: }

2138: /*@
2139:    TSSetTime - Allows one to reset the time.

2141:    Logically Collective on TS

2143:    Input Parameters:
2144: +  ts - the TS context obtained from TSCreate()
2145: -  time - the time

2147:    Level: intermediate

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

2151: .keywords: TS, set, time
2152: @*/
2153: PetscErrorCode  TSSetTime(TS ts, PetscReal t)
2154: {
2158:   ts->ptime = t;
2159:   return(0);
2160: }

2164: /*@C
2165:    TSSetOptionsPrefix - Sets the prefix used for searching for all
2166:    TS options in the database.

2168:    Logically Collective on TS

2170:    Input Parameter:
2171: +  ts     - The TS context
2172: -  prefix - The prefix to prepend to all option names

2174:    Notes:
2175:    A hyphen (-) must NOT be given at the beginning of the prefix name.
2176:    The first character of all runtime options is AUTOMATICALLY the
2177:    hyphen.

2179:    Level: advanced

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

2183: .seealso: TSSetFromOptions()

2185: @*/
2186: PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
2187: {
2189:   SNES           snes;

2193:   PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);
2194:   TSGetSNES(ts,&snes);
2195:   SNESSetOptionsPrefix(snes,prefix);
2196:   return(0);
2197: }


2202: /*@C
2203:    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
2204:    TS options in the database.

2206:    Logically Collective on TS

2208:    Input Parameter:
2209: +  ts     - The TS context
2210: -  prefix - The prefix to prepend to all option names

2212:    Notes:
2213:    A hyphen (-) must NOT be given at the beginning of the prefix name.
2214:    The first character of all runtime options is AUTOMATICALLY the
2215:    hyphen.

2217:    Level: advanced

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

2221: .seealso: TSGetOptionsPrefix()

2223: @*/
2224: PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
2225: {
2227:   SNES           snes;

2231:   PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);
2232:   TSGetSNES(ts,&snes);
2233:   SNESAppendOptionsPrefix(snes,prefix);
2234:   return(0);
2235: }

2239: /*@C
2240:    TSGetOptionsPrefix - Sets the prefix used for searching for all
2241:    TS options in the database.

2243:    Not Collective

2245:    Input Parameter:
2246: .  ts - The TS context

2248:    Output Parameter:
2249: .  prefix - A pointer to the prefix string used

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

2254:    Level: intermediate

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

2258: .seealso: TSAppendOptionsPrefix()
2259: @*/
2260: PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
2261: {

2267:   PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);
2268:   return(0);
2269: }

2273: /*@C
2274:    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.

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

2278:    Input Parameter:
2279: .  ts  - The TS context obtained from TSCreate()

2281:    Output Parameters:
2282: +  J   - The Jacobian J of F, where U_t = F(U,t)
2283: .  M   - The preconditioner matrix, usually the same as J
2284: .  func - Function to compute the Jacobian of the RHS
2285: -  ctx - User-defined context for Jacobian evaluation routine

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

2289:    Level: intermediate

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

2293: .keywords: TS, timestep, get, matrix, Jacobian
2294: @*/
2295: PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx)
2296: {
2298:   SNES           snes;

2301:   TSGetSNES(ts,&snes);
2302:   SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);
2303:   if (func) *func = ts->userops->rhsjacobian;
2304:   if (ctx) *ctx = ts->jacP;
2305:   return(0);
2306: }

2310: /*@C
2311:    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.

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

2315:    Input Parameter:
2316: .  ts  - The TS context obtained from TSCreate()

2318:    Output Parameters:
2319: +  A   - The Jacobian of F(t,U,U_t)
2320: .  B   - The preconditioner matrix, often the same as A
2321: .  f   - The function to compute the matrices
2322: - ctx - User-defined context for Jacobian evaluation routine

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

2326:    Level: advanced

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

2330: .keywords: TS, timestep, get, matrix, Jacobian
2331: @*/
2332: PetscErrorCode  TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
2333: {
2335:   SNES           snes;

2338:   TSGetSNES(ts,&snes);
2339:   SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);
2340:   if (f) *f = ts->userops->ijacobian;
2341:   if (ctx) *ctx = ts->jacP;
2342:   return(0);
2343: }

2345: typedef struct {
2346:   PetscViewer viewer;
2347:   Vec         initialsolution;
2348:   PetscBool   showinitial;
2349: } TSMonitorSolutionCtx;

2353: /*@C
2354:    TSMonitorSolution - Monitors progress of the TS solvers by calling 
2355:    VecView() for the solution at each timestep

2357:    Collective on TS

2359:    Input Parameters:
2360: +  ts - the TS context
2361: .  step - current time-step
2362: .  ptime - current time
2363: -  dummy - either a viewer or PETSC_NULL

2365:    Level: intermediate

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

2369: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
2370: @*/
2371: PetscErrorCode  TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
2372: {
2373:   PetscErrorCode       ierr;
2374:   TSMonitorSolutionCtx *ictx = (TSMonitorSolutionCtx*)dummy;

2377:   if (!step && ictx->showinitial) {
2378:     if (!ictx->initialsolution) {
2379:       VecDuplicate(x,&ictx->initialsolution);
2380:     }
2381:     VecCopy(x,ictx->initialsolution);
2382:   }
2383:   if (ictx->showinitial) {
2384:     PetscReal pause;
2385:     PetscViewerDrawGetPause(ictx->viewer,&pause);
2386:     PetscViewerDrawSetPause(ictx->viewer,0.0);
2387:     VecView(ictx->initialsolution,ictx->viewer);
2388:     PetscViewerDrawSetPause(ictx->viewer,pause);
2389:     PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);
2390:   }
2391:   VecView(x,ictx->viewer);
2392:   if (ictx->showinitial) {
2393:     PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);
2394:   }
2395:   return(0);
2396: }


2401: /*@C
2402:    TSMonitorSolutionDestroy - Destroys the monitor context for TSMonitorSolution

2404:    Collective on TS

2406:    Input Parameters:
2407: .    ctx - the monitor context

2409:    Level: intermediate

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

2413: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
2414: @*/
2415: PetscErrorCode  TSMonitorSolutionDestroy(void **ctx)
2416: {
2417:   PetscErrorCode       ierr;
2418:   TSMonitorSolutionCtx *ictx = *(TSMonitorSolutionCtx**)ctx;
2419: 
2421:   PetscViewerDestroy(&ictx->viewer);
2422:   VecDestroy(&ictx->initialsolution);
2423:   PetscFree(ictx);
2424:   return(0);
2425: }

2429: /*@C
2430:    TSMonitorSolutionCreate - Creates the monitor context for TSMonitorSolution

2432:    Collective on TS

2434:    Input Parameter:
2435: .    ts - time-step context

2437:    Output Patameter:
2438: .    ctx - the monitor context

2440:    Level: intermediate

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

2444: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
2445: @*/
2446: PetscErrorCode  TSMonitorSolutionCreate(TS ts,PetscViewer viewer,void **ctx)
2447: {
2448:   PetscErrorCode       ierr;
2449:   TSMonitorSolutionCtx *ictx;
2450: 
2452:   PetscNew(TSMonitorSolutionCtx,&ictx);
2453:   *ctx = (void*)ictx;
2454:   if (!viewer) {
2455:     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
2456:   }
2457:   PetscObjectReference((PetscObject)viewer);
2458:   ictx->viewer      = viewer;
2459:   ictx->showinitial = PETSC_FALSE;
2460:   PetscOptionsGetBool(((PetscObject)ts)->prefix,"-ts_monitor_solution_initial",&ictx->showinitial,PETSC_NULL);
2461:   return(0);
2462: }

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

2469:    Logically Collective on TS and DM

2471:    Input Parameters:
2472: +  ts - the preconditioner context
2473: -  dm - the dm

2475:    Level: intermediate


2478: .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
2479: @*/
2480: PetscErrorCode  TSSetDM(TS ts,DM dm)
2481: {
2483:   SNES           snes;

2487:   PetscObjectReference((PetscObject)dm);
2488:   DMDestroy(&ts->dm);
2489:   ts->dm = dm;
2490:   TSGetSNES(ts,&snes);
2491:   SNESSetDM(snes,dm);
2492:   return(0);
2493: }

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

2500:    Not Collective

2502:    Input Parameter:
2503: . ts - the preconditioner context

2505:    Output Parameter:
2506: .  dm - the dm

2508:    Level: intermediate


2511: .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
2512: @*/
2513: PetscErrorCode  TSGetDM(TS ts,DM *dm)
2514: {

2519:   if (!ts->dm) {
2520:     DMShellCreate(((PetscObject)ts)->comm,&ts->dm);
2521:     if (ts->snes) {SNESSetDM(ts->snes,ts->dm);}
2522:   }
2523:   *dm = ts->dm;
2524:   return(0);
2525: }

2529: /*@
2530:    SNESTSFormFunction - Function to evaluate nonlinear residual

2532:    Logically Collective on SNES

2534:    Input Parameter:
2535: + snes - nonlinear solver
2536: . X - the current state at which to evaluate the residual
2537: - ctx - user context, must be a TS

2539:    Output Parameter:
2540: . F - the nonlinear residual

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

2546:    Level: advanced

2548: .seealso: SNESSetFunction(), MatFDColoringSetFunction()
2549: @*/
2550: PetscErrorCode  SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
2551: {
2552:   TS ts = (TS)ctx;

2560:   (ts->ops->snesfunction)(snes,X,F,ts);
2561:   return(0);
2562: }

2566: /*@
2567:    SNESTSFormJacobian - Function to evaluate the Jacobian

2569:    Collective on SNES

2571:    Input Parameter:
2572: + snes - nonlinear solver
2573: . X - the current state at which to evaluate the residual
2574: - ctx - user context, must be a TS

2576:    Output Parameter:
2577: + A - the Jacobian
2578: . B - the preconditioning matrix (may be the same as A)
2579: - flag - indicates any structure change in the matrix

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

2584:    Level: developer

2586: .seealso: SNESSetJacobian()
2587: @*/
2588: PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
2589: {
2590:   TS ts = (TS)ctx;

2602:   (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);
2603:   return(0);
2604: }

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

2611:    Collective on TS

2613:    Input Arguments:
2614: +  ts - time stepping context
2615: .  t - time at which to evaluate
2616: .  X - state at which to evaluate
2617: -  ctx - context

2619:    Output Arguments:
2620: .  F - right hand side

2622:    Level: intermediate

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

2628: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
2629: @*/
2630: PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
2631: {
2633:   Mat Arhs,Brhs;
2634:   MatStructure flg2;

2637:   TSGetRHSMats_Private(ts,&Arhs,&Brhs);
2638:   TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);
2639:   MatMult(Arhs,X,F);
2640:   return(0);
2641: }

2645: /*@C
2646:    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.

2648:    Collective on TS

2650:    Input Arguments:
2651: +  ts - time stepping context
2652: .  t - time at which to evaluate
2653: .  X - state at which to evaluate
2654: -  ctx - context

2656:    Output Arguments:
2657: +  A - pointer to operator
2658: .  B - pointer to preconditioning matrix
2659: -  flg - matrix structure flag

2661:    Level: intermediate

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

2666: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
2667: @*/
2668: PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg,void *ctx)
2669: {

2672:   *flg = SAME_PRECONDITIONER;
2673:   return(0);
2674: }

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

2681:    Collective on TS

2683:    Input Arguments:
2684: +  ts - time stepping context
2685: .  t - time at which to evaluate
2686: .  X - state at which to evaluate
2687: .  Xdot - time derivative of state vector
2688: -  ctx - context

2690:    Output Arguments:
2691: .  F - left hand side

2693:    Level: intermediate

2695:    Notes:
2696:    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
2697:    user is required to write their own TSComputeIFunction.
2698:    This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems.
2699:    The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian().

2701: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant()
2702: @*/
2703: PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
2704: {
2706:   Mat A,B;
2707:   MatStructure flg2;

2710:   TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);
2711:   TSComputeIJacobian(ts,t,X,Xdot,1.0,&A,&B,&flg2,PETSC_TRUE);
2712:   MatMult(A,Xdot,F);
2713:   return(0);
2714: }

2718: /*@C
2719:    TSComputeIJacobianConstant - Reuses a Jacobian that is time-independent.

2721:    Collective on TS

2723:    Input Arguments:
2724: +  ts - time stepping context
2725: .  t - time at which to evaluate
2726: .  X - state at which to evaluate
2727: .  Xdot - time derivative of state vector
2728: .  shift - shift to apply
2729: -  ctx - context

2731:    Output Arguments:
2732: +  A - pointer to operator
2733: .  B - pointer to preconditioning matrix
2734: -  flg - matrix structure flag

2736:    Level: intermediate

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

2741: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
2742: @*/
2743: PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx)
2744: {

2747:   *flg = SAME_PRECONDITIONER;
2748:   return(0);
2749: }


2754: /*@
2755:    TSGetConvergedReason - Gets the reason the TS iteration was stopped.

2757:    Not Collective

2759:    Input Parameter:
2760: .  ts - the TS context

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

2766:    Level: intermediate

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

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

2773: .seealso: TSSetConvergenceTest(), TSConvergedReason
2774: @*/
2775: PetscErrorCode  TSGetConvergedReason(TS ts,TSConvergedReason *reason)
2776: {
2780:   *reason = ts->reason;
2781:   return(0);
2782: }

2786: /*@
2787:    TSGetSNESIterations - Gets the total number of nonlinear iterations
2788:    used by the time integrator.

2790:    Not Collective

2792:    Input Parameter:
2793: .  ts - TS context

2795:    Output Parameter:
2796: .  nits - number of nonlinear iterations

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

2801:    Level: intermediate

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

2805: .seealso:  TSGetKSPIterations()
2806: @*/
2807: PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits)
2808: {
2812:   *nits = ts->snes_its;
2813:   return(0);
2814: }

2818: /*@
2819:    TSGetKSPIterations - Gets the total number of linear iterations
2820:    used by the time integrator.

2822:    Not Collective

2824:    Input Parameter:
2825: .  ts - TS context

2827:    Output Parameter:
2828: .  lits - number of linear iterations

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

2833:    Level: intermediate

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

2837: .seealso:  TSGetSNESIterations(), SNESGetKSPIterations()
2838: @*/
2839: PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits)
2840: {
2844:   *lits = ts->ksp_its;
2845:   return(0);
2846: }

2850: /*@
2851:    TSGetStepRejections - Gets the total number of rejected steps.

2853:    Not Collective

2855:    Input Parameter:
2856: .  ts - TS context

2858:    Output Parameter:
2859: .  rejects - number of steps rejected

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

2864:    Level: intermediate

2866: .keywords: TS, get, number

2868: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails()
2869: @*/
2870: PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects)
2871: {
2875:   *rejects = ts->reject;
2876:   return(0);
2877: }

2881: /*@
2882:    TSGetSNESFailures - Gets the total number of failed SNES solves

2884:    Not Collective

2886:    Input Parameter:
2887: .  ts - TS context

2889:    Output Parameter:
2890: .  fails - number of failed nonlinear solves

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

2895:    Level: intermediate

2897: .keywords: TS, get, number

2899: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures()
2900: @*/
2901: PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails)
2902: {
2906:   *fails = ts->num_snes_failures;
2907:   return(0);
2908: }

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

2915:    Not Collective

2917:    Input Parameter:
2918: +  ts - TS context
2919: -  rejects - maximum number of rejected steps, pass -1 for unlimited

2921:    Notes:
2922:    The counter is reset to zero for each step

2924:    Options Database Key:
2925:  .  -ts_max_reject - Maximum number of step rejections before a step fails

2927:    Level: intermediate

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

2931: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
2932: @*/
2933: PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects)
2934: {
2937:   ts->max_reject = rejects;
2938:   return(0);
2939: }

2943: /*@
2944:    TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves

2946:    Not Collective

2948:    Input Parameter:
2949: +  ts - TS context
2950: -  fails - maximum number of failed nonlinear solves, pass -1 for unlimited

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

2955:    Options Database Key:
2956:  .  -ts_max_snes_failures - Maximum number of nonlinear solve failures

2958:    Level: intermediate

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

2962: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason()
2963: @*/
2964: PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails)
2965: {
2968:   ts->max_snes_failures = fails;
2969:   return(0);
2970: }

2974: /*@
2975:    TSSetErrorIfStepFails - Error if no step succeeds

2977:    Not Collective

2979:    Input Parameter:
2980: +  ts - TS context
2981: -  err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure

2983:    Options Database Key:
2984:  .  -ts_error_if_step_fails - Error if no step succeeds

2986:    Level: intermediate

2988: .keywords: TS, set, error

2990: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
2991: @*/
2992: PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err)
2993: {
2996:   ts->errorifstepfailed = err;
2997:   return(0);
2998: }

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

3005:    Collective on TS

3007:    Input Parameters:
3008: +  ts - the TS context
3009: .  step - current time-step
3010: .  ptime - current time
3011: .  x - current state
3012: -  viewer - binary viewer

3014:    Level: intermediate

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

3018: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3019: @*/
3020: PetscErrorCode  TSMonitorSolutionBinary(TS ts,PetscInt step,PetscReal ptime,Vec x,void *viewer)
3021: {
3022:   PetscErrorCode       ierr;
3023:   PetscViewer          v = (PetscViewer)viewer;

3026:   VecView(x,v);
3027:   return(0);
3028: }

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

3035:    Collective on TS

3037:    Input Parameters:
3038: +  ts - the TS context
3039: .  step - current time-step
3040: .  ptime - current time
3041: .  x - current state
3042: -  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)

3044:    Level: intermediate

3046:    Notes:
3047:    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.
3048:    These are named according to the file name template.

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

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

3054: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3055: @*/
3056: PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec x,void *filenametemplate)
3057: {
3059:   char           filename[PETSC_MAX_PATH_LEN];
3060:   PetscViewer    viewer;

3063:   PetscSNPrintf(filename,sizeof filename,(const char*)filenametemplate,step);
3064:   PetscViewerVTKOpen(((PetscObject)ts)->comm,filename,FILE_MODE_WRITE,&viewer);
3065:   VecView(x,viewer);
3066:   PetscViewerDestroy(&viewer);
3067:   return(0);
3068: }

3072: /*@C
3073:    TSMonitorSolutionVTKDestroy - Destroy context for monitoring

3075:    Collective on TS

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

3080:    Level: intermediate

3082:    Note:
3083:    This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK().

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

3087: .seealso: TSMonitorSet(), TSMonitorSolutionVTK()
3088: @*/
3089: PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
3090: {

3094:   PetscFree(*(char**)filenametemplate);
3095:   return(0);
3096: }

3100: /*@
3101:    TSGetAdapt - Get the adaptive controller context for the current method

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

3105:    Input Arguments:
3106: .  ts - time stepping context

3108:    Output Arguments:
3109: .  adapt - adaptive controller

3111:    Level: intermediate

3113: .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose()
3114: @*/
3115: PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt)
3116: {

3122:   if (!ts->adapt) {
3123:     TSAdaptCreate(((PetscObject)ts)->comm,&ts->adapt);
3124:     PetscLogObjectParent(ts,ts->adapt);
3125:     PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);
3126:   }
3127:   *adapt = ts->adapt;
3128:   return(0);
3129: }

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

3136:    Logically Collective

3138:    Input Arguments:
3139: +  ts - time integration context
3140: .  atol - scalar absolute tolerances, PETSC_DECIDE to leave current value
3141: .  vatol - vector of absolute tolerances or PETSC_NULL, used in preference to atol if present
3142: .  rtol - scalar relative tolerances, PETSC_DECIDE to leave current value
3143: -  vrtol - vector of relative tolerances or PETSC_NULL, used in preference to atol if present

3145:    Level: beginner

3147: .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances()
3148: @*/
3149: PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol)
3150: {

3154:   if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol;
3155:   if (vatol) {
3156:     PetscObjectReference((PetscObject)vatol);
3157:     VecDestroy(&ts->vatol);
3158:     ts->vatol = vatol;
3159:   }
3160:   if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol;
3161:   if (vrtol) {
3162:     PetscObjectReference((PetscObject)vrtol);
3163:     VecDestroy(&ts->vrtol);
3164:     ts->vrtol = vrtol;
3165:   }
3166:   return(0);
3167: }

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

3174:    Logically Collective

3176:    Input Arguments:
3177: .  ts - time integration context

3179:    Output Arguments:
3180: +  atol - scalar absolute tolerances, PETSC_NULL to ignore
3181: .  vatol - vector of absolute tolerances, PETSC_NULL to ignore
3182: .  rtol - scalar relative tolerances, PETSC_NULL to ignore
3183: -  vrtol - vector of relative tolerances, PETSC_NULL to ignore

3185:    Level: beginner

3187: .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances()
3188: @*/
3189: PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol)
3190: {

3193:   if (atol)  *atol  = ts->atol;
3194:   if (vatol) *vatol = ts->vatol;
3195:   if (rtol)  *rtol  = ts->rtol;
3196:   if (vrtol) *vrtol = ts->vrtol;
3197:   return(0);
3198: }

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

3205:    Collective on TS

3207:    Input Arguments:
3208: +  ts - time stepping context
3209: -  Y - state vector to be compared to ts->vec_sol

3211:    Output Arguments:
3212: .  norm - weighted norm, a value of 1.0 is considered small

3214:    Level: developer

3216: .seealso: TSSetTolerances()
3217: @*/
3218: PetscErrorCode TSErrorNormWRMS(TS ts,Vec Y,PetscReal *norm)
3219: {
3221:   PetscInt i,n,N;
3222:   const PetscScalar *x,*y;
3223:   Vec X;
3224:   PetscReal sum,gsum;

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

3234:   VecGetSize(X,&N);
3235:   VecGetLocalSize(X,&n);
3236:   VecGetArrayRead(X,&x);
3237:   VecGetArrayRead(Y,&y);
3238:   sum = 0.;
3239:   if (ts->vatol && ts->vrtol) {
3240:     const PetscScalar *atol,*rtol;
3241:     VecGetArrayRead(ts->vatol,&atol);
3242:     VecGetArrayRead(ts->vrtol,&rtol);
3243:     for (i=0; i<n; i++) {
3244:       PetscReal tol = PetscRealPart(atol[i]) + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3245:       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3246:     }
3247:     VecRestoreArrayRead(ts->vatol,&atol);
3248:     VecRestoreArrayRead(ts->vrtol,&rtol);
3249:   } else if (ts->vatol) {       /* vector atol, scalar rtol */
3250:     const PetscScalar *atol;
3251:     VecGetArrayRead(ts->vatol,&atol);
3252:     for (i=0; i<n; i++) {
3253:       PetscReal tol = PetscRealPart(atol[i]) + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3254:       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3255:     }
3256:     VecRestoreArrayRead(ts->vatol,&atol);
3257:   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
3258:     const PetscScalar *rtol;
3259:     VecGetArrayRead(ts->vrtol,&rtol);
3260:     for (i=0; i<n; i++) {
3261:       PetscReal tol = ts->atol + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3262:       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3263:     }
3264:     VecRestoreArrayRead(ts->vrtol,&rtol);
3265:   } else {                      /* scalar atol, scalar rtol */
3266:     for (i=0; i<n; i++) {
3267:       PetscReal tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3268:       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3269:     }
3270:   }
3271:   VecRestoreArrayRead(X,&x);
3272:   VecRestoreArrayRead(Y,&y);

3274:   MPI_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,((PetscObject)ts)->comm);
3275:   *norm = PetscSqrtReal(gsum / N);
3276:   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
3277:   return(0);
3278: }

3282: /*@
3283:    TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler

3285:    Logically Collective on TS

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

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

3294:    Level: intermediate

3296: .seealso: TSGetCFLTime(), TSADAPTCFL
3297: @*/
3298: PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime)
3299: {

3303:   ts->cfltime_local = cfltime;
3304:   ts->cfltime = -1.;
3305:   return(0);
3306: }

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

3313:    Collective on TS

3315:    Input Arguments:
3316: .  ts - time stepping context

3318:    Output Arguments:
3319: .  cfltime - maximum stable time step for forward Euler

3321:    Level: advanced

3323: .seealso: TSSetCFLTimeLocal()
3324: @*/
3325: PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime)
3326: {

3330:   if (ts->cfltime < 0) {
3331:     MPI_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,((PetscObject)ts)->comm);
3332:   }
3333:   *cfltime = ts->cfltime;
3334:   return(0);
3335: }

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

3342:    Input Parameters:
3343: .  ts   - the TS context.
3344: .  xl   - lower bound.
3345: .  xu   - upper bound.

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

3351:    Level: advanced

3353: @*/
3354: PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
3355: {
3357:   SNES           snes;

3360:   TSGetSNES(ts,&snes);
3361:   SNESVISetVariableBounds(snes,xl,xu);
3362:   return(0);
3363: }

3365: #if defined(PETSC_HAVE_MATLAB_ENGINE)
3366: #include <mex.h>

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

3372: /*
3373:    TSComputeFunction_Matlab - Calls the function that has been set with
3374:                          TSSetFunctionMatlab().

3376:    Collective on TS

3378:    Input Parameters:
3379: +  snes - the TS context
3380: -  x - input vector

3382:    Output Parameter:
3383: .  y - function vector, as set by TSSetFunction()

3385:    Notes:
3386:    TSComputeFunction() is typically used within nonlinear solvers
3387:    implementations, so most users would not generally call this routine
3388:    themselves.

3390:    Level: developer

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

3394: .seealso: TSSetFunction(), TSGetFunction()
3395: */
3396: PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx)
3397: {
3398:   PetscErrorCode   ierr;
3399:   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3400:   int              nlhs = 1,nrhs = 7;
3401:   mxArray          *plhs[1],*prhs[7];
3402:   long long int    lx = 0,lxdot = 0,ly = 0,ls = 0;


3412:   PetscMemcpy(&ls,&snes,sizeof(snes));
3413:   PetscMemcpy(&lx,&x,sizeof(x));
3414:   PetscMemcpy(&lxdot,&xdot,sizeof(xdot));
3415:   PetscMemcpy(&ly,&y,sizeof(x));
3416:   prhs[0] =  mxCreateDoubleScalar((double)ls);
3417:   prhs[1] =  mxCreateDoubleScalar(time);
3418:   prhs[2] =  mxCreateDoubleScalar((double)lx);
3419:   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
3420:   prhs[4] =  mxCreateDoubleScalar((double)ly);
3421:   prhs[5] =  mxCreateString(sctx->funcname);
3422:   prhs[6] =  sctx->ctx;
3423:    mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");
3424:    mxGetScalar(plhs[0]);
3425:   mxDestroyArray(prhs[0]);
3426:   mxDestroyArray(prhs[1]);
3427:   mxDestroyArray(prhs[2]);
3428:   mxDestroyArray(prhs[3]);
3429:   mxDestroyArray(prhs[4]);
3430:   mxDestroyArray(prhs[5]);
3431:   mxDestroyArray(plhs[0]);
3432:   return(0);
3433: }


3438: /*
3439:    TSSetFunctionMatlab - Sets the function evaluation routine and function
3440:    vector for use by the TS routines in solving ODEs
3441:    equations from MATLAB. Here the function is a string containing the name of a MATLAB function

3443:    Logically Collective on TS

3445:    Input Parameters:
3446: +  ts - the TS context
3447: -  func - function evaluation routine

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

3452:    Level: beginner

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

3456: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3457: */
3458: PetscErrorCode  TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx)
3459: {
3460:   PetscErrorCode  ierr;
3461:   TSMatlabContext *sctx;

3464:   /* currently sctx is memory bleed */
3465:   PetscMalloc(sizeof(TSMatlabContext),&sctx);
3466:   PetscStrallocpy(func,&sctx->funcname);
3467:   /*
3468:      This should work, but it doesn't
3469:   sctx->ctx = ctx;
3470:   mexMakeArrayPersistent(sctx->ctx);
3471:   */
3472:   sctx->ctx = mxDuplicateArray(ctx);
3473:   TSSetIFunction(ts,PETSC_NULL,TSComputeFunction_Matlab,sctx);
3474:   return(0);
3475: }

3479: /*
3480:    TSComputeJacobian_Matlab - Calls the function that has been set with
3481:                          TSSetJacobianMatlab().

3483:    Collective on TS

3485:    Input Parameters:
3486: +  ts - the TS context
3487: .  x - input vector
3488: .  A, B - the matrices
3489: -  ctx - user context

3491:    Output Parameter:
3492: .  flag - structure of the matrix

3494:    Level: developer

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

3498: .seealso: TSSetFunction(), TSGetFunction()
3499: @*/
3500: PetscErrorCode  TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx)
3501: {
3502:   PetscErrorCode  ierr;
3503:   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3504:   int             nlhs = 2,nrhs = 9;
3505:   mxArray         *plhs[2],*prhs[9];
3506:   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;


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

3514:   PetscMemcpy(&ls,&ts,sizeof(ts));
3515:   PetscMemcpy(&lx,&x,sizeof(x));
3516:   PetscMemcpy(&lxdot,&xdot,sizeof(x));
3517:   PetscMemcpy(&lA,A,sizeof(x));
3518:   PetscMemcpy(&lB,B,sizeof(x));
3519:   prhs[0] =  mxCreateDoubleScalar((double)ls);
3520:   prhs[1] =  mxCreateDoubleScalar((double)time);
3521:   prhs[2] =  mxCreateDoubleScalar((double)lx);
3522:   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
3523:   prhs[4] =  mxCreateDoubleScalar((double)shift);
3524:   prhs[5] =  mxCreateDoubleScalar((double)lA);
3525:   prhs[6] =  mxCreateDoubleScalar((double)lB);
3526:   prhs[7] =  mxCreateString(sctx->funcname);
3527:   prhs[8] =  sctx->ctx;
3528:    mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");
3529:    mxGetScalar(plhs[0]);
3530:   *flag   =  (MatStructure) mxGetScalar(plhs[1]);
3531:   mxDestroyArray(prhs[0]);
3532:   mxDestroyArray(prhs[1]);
3533:   mxDestroyArray(prhs[2]);
3534:   mxDestroyArray(prhs[3]);
3535:   mxDestroyArray(prhs[4]);
3536:   mxDestroyArray(prhs[5]);
3537:   mxDestroyArray(prhs[6]);
3538:   mxDestroyArray(prhs[7]);
3539:   mxDestroyArray(plhs[0]);
3540:   mxDestroyArray(plhs[1]);
3541:   return(0);
3542: }


3547: /*
3548:    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
3549:    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

3551:    Logically Collective on TS

3553:    Input Parameters:
3554: +  ts - the TS context
3555: .  A,B - Jacobian matrices
3556: .  func - function evaluation routine
3557: -  ctx - user context

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


3563:    Level: developer

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

3567: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3568: */
3569: PetscErrorCode  TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx)
3570: {
3571:   PetscErrorCode    ierr;
3572:   TSMatlabContext *sctx;

3575:   /* currently sctx is memory bleed */
3576:   PetscMalloc(sizeof(TSMatlabContext),&sctx);
3577:   PetscStrallocpy(func,&sctx->funcname);
3578:   /*
3579:      This should work, but it doesn't
3580:   sctx->ctx = ctx;
3581:   mexMakeArrayPersistent(sctx->ctx);
3582:   */
3583:   sctx->ctx = mxDuplicateArray(ctx);
3584:   TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);
3585:   return(0);
3586: }

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

3593:    Collective on TS

3595: .seealso: TSSetFunction(), TSGetFunction()
3596: @*/
3597: PetscErrorCode  TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec x, void *ctx)
3598: {
3599:   PetscErrorCode  ierr;
3600:   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3601:   int             nlhs = 1,nrhs = 6;
3602:   mxArray         *plhs[1],*prhs[6];
3603:   long long int   lx = 0,ls = 0;


3609:   PetscMemcpy(&ls,&ts,sizeof(ts));
3610:   PetscMemcpy(&lx,&x,sizeof(x));
3611:   prhs[0] =  mxCreateDoubleScalar((double)ls);
3612:   prhs[1] =  mxCreateDoubleScalar((double)it);
3613:   prhs[2] =  mxCreateDoubleScalar((double)time);
3614:   prhs[3] =  mxCreateDoubleScalar((double)lx);
3615:   prhs[4] =  mxCreateString(sctx->funcname);
3616:   prhs[5] =  sctx->ctx;
3617:    mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");
3618:    mxGetScalar(plhs[0]);
3619:   mxDestroyArray(prhs[0]);
3620:   mxDestroyArray(prhs[1]);
3621:   mxDestroyArray(prhs[2]);
3622:   mxDestroyArray(prhs[3]);
3623:   mxDestroyArray(prhs[4]);
3624:   mxDestroyArray(plhs[0]);
3625:   return(0);
3626: }


3631: /*
3632:    TSMonitorSetMatlab - Sets the monitor function from Matlab

3634:    Level: developer

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

3638: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3639: */
3640: PetscErrorCode  TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx)
3641: {
3642:   PetscErrorCode    ierr;
3643:   TSMatlabContext *sctx;

3646:   /* currently sctx is memory bleed */
3647:   PetscMalloc(sizeof(TSMatlabContext),&sctx);
3648:   PetscStrallocpy(func,&sctx->funcname);
3649:   /*
3650:      This should work, but it doesn't
3651:   sctx->ctx = ctx;
3652:   mexMakeArrayPersistent(sctx->ctx);
3653:   */
3654:   sctx->ctx = mxDuplicateArray(ctx);
3655:   TSMonitorSet(ts,TSMonitor_Matlab,sctx,PETSC_NULL);
3656:   return(0);
3657: }
3658: #endif