TSIRKTableauCreate#

create the tableau for TSIRK and provide the entries

Synopsis#

#include "petscts.h"   
PetscErrorCode TSIRKTableauCreate(TS ts, PetscInt nstages, const PetscReal *A, const PetscReal *b, const PetscReal *c, const PetscReal *binterp, const PetscScalar *A_inv, const PetscScalar *A_inv_rowsum, const PetscScalar *I_s)

Not Collective

Input Parameters#

  • ts - timestepping context

  • nstages - number of stages, this is the dimension of the matrices below

  • A - stage coefficients (dimension nstages*nstages, row-major)

  • b - step completion table (dimension nstages)

  • c - abscissa (dimension nstages)

  • binterp - coefficients of the interpolation formula (dimension nstages)

  • A_inv - inverse of A (dimension nstages*nstages, row-major)

  • A_inv_rowsum - row sum of the inverse of A (dimension nstages)

  • I_s - identity matrix (dimension nstages*nstages)

See Also#

TS: Scalable ODE and DAE Solvers, TSIRK, TSIRKRegister()

Level#

advanced

Location#

src/ts/impls/implicit/irk/irk.c


Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages