xref: /petsc/src/ts/impls/implicit/radau5/radau5.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1d249a78fSBarry Smith /*
2d249a78fSBarry Smith     Provides a PETSc interface to RADAU5 solver.
3d249a78fSBarry Smith 
4d249a78fSBarry Smith */
5d249a78fSBarry Smith #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
6d249a78fSBarry Smith 
7d249a78fSBarry Smith typedef struct {
86ffeb343SBarry Smith   Vec work, workf;
9d249a78fSBarry Smith } TS_Radau5;
10d249a78fSBarry Smith 
11*d71ae5a4SJacob Faibussowitsch void FVPOL(int *N, double *X, double *Y, double *F, double *RPAR, void *IPAR)
12*d71ae5a4SJacob Faibussowitsch {
136ffeb343SBarry Smith   TS          ts    = (TS)IPAR;
146ffeb343SBarry Smith   TS_Radau5  *cvode = (TS_Radau5 *)ts->data;
156ffeb343SBarry Smith   DM          dm;
166ffeb343SBarry Smith   DMTS        tsdm;
176ffeb343SBarry Smith   TSIFunction ifunction;
186ffeb343SBarry Smith 
199566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, VecPlaceArray(cvode->work, Y));
209566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, VecPlaceArray(cvode->workf, F));
216ffeb343SBarry Smith 
226ffeb343SBarry Smith   /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */
239566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, TSGetDM(ts, &dm));
249566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, DMGetDMTS(dm, &tsdm));
259566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, DMTSGetIFunction(dm, &ifunction, NULL));
266ffeb343SBarry Smith   if (!ifunction) {
279566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_SELF, TSComputeRHSFunction(ts, *X, cvode->work, cvode->workf));
28e68ccae8SBarry Smith   } else { /* If rhsfunction is also set, this computes both parts and scale them to the right hand side */
296ffeb343SBarry Smith     Vec yydot;
306ffeb343SBarry Smith 
319566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_SELF, VecDuplicate(cvode->work, &yydot));
329566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_SELF, VecZeroEntries(yydot));
339566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_SELF, TSComputeIFunction(ts, *X, cvode->work, yydot, cvode->workf, PETSC_FALSE));
349566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_SELF, VecScale(cvode->workf, -1.));
359566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_SELF, VecDestroy(&yydot));
366ffeb343SBarry Smith   }
376ffeb343SBarry Smith 
389566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, VecResetArray(cvode->work));
399566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, VecResetArray(cvode->workf));
406ffeb343SBarry Smith }
416ffeb343SBarry Smith 
42*d71ae5a4SJacob Faibussowitsch void JVPOL(PetscInt *N, PetscScalar *X, PetscScalar *Y, PetscScalar *DFY, int *LDFY, PetscScalar *RPAR, void *IPAR)
43*d71ae5a4SJacob Faibussowitsch {
448acf3572SBarry Smith   TS         ts    = (TS)IPAR;
458acf3572SBarry Smith   TS_Radau5 *cvode = (TS_Radau5 *)ts->data;
468acf3572SBarry Smith   Vec        yydot;
478acf3572SBarry Smith   Mat        mat;
488acf3572SBarry Smith   PetscInt   n;
49d249a78fSBarry Smith 
509566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, VecPlaceArray(cvode->work, Y));
519566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, VecDuplicate(cvode->work, &yydot));
529566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, VecGetSize(yydot, &n));
539566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, MatCreateSeqDense(PETSC_COMM_SELF, n, n, DFY, &mat));
549566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, VecZeroEntries(yydot));
559566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, TSComputeIJacobian(ts, *X, cvode->work, yydot, 0, mat, mat, PETSC_FALSE));
569566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, MatScale(mat, -1.0));
579566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, MatDestroy(&mat));
589566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, VecDestroy(&yydot));
599566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, VecResetArray(cvode->work));
608acf3572SBarry Smith }
618acf3572SBarry Smith 
62*d71ae5a4SJacob Faibussowitsch void SOLOUT(int *NR, double *XOLD, double *X, double *Y, double *CONT, double *LRC, int *N, double *RPAR, void *IPAR, int *IRTRN)
63*d71ae5a4SJacob Faibussowitsch {
64d249a78fSBarry Smith   TS         ts    = (TS)IPAR;
65d249a78fSBarry Smith   TS_Radau5 *cvode = (TS_Radau5 *)ts->data;
66d249a78fSBarry Smith 
679566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, VecPlaceArray(cvode->work, Y));
68d249a78fSBarry Smith   ts->time_step = *X - *XOLD;
699566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, TSMonitor(ts, *NR - 1, *X, cvode->work));
709566063dSJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, VecResetArray(cvode->work));
71d249a78fSBarry Smith }
72d249a78fSBarry Smith 
73d249a78fSBarry Smith void radau5_(int *, void *, double *, double *, double *, double *, double *, double *, int *, void *, int *, int *, int *, void *, int *, int *, int *, void *, int *, double *, int *, int *, int *, double *, void *, int *);
74d249a78fSBarry Smith 
75*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSolve_Radau5(TS ts)
76*d71ae5a4SJacob Faibussowitsch {
77d249a78fSBarry Smith   TS_Radau5   *cvode = (TS_Radau5 *)ts->data;
78d249a78fSBarry Smith   PetscScalar *Y, *WORK, X, XEND, RTOL, ATOL, H, RPAR;
79d249a78fSBarry Smith   PetscInt     ND, *IWORK, LWORK, LIWORK, MUJAC, MLMAS, MUMAS, IDID, ITOL;
80d249a78fSBarry Smith   int          IJAC, MLJAC, IMAS, IOUT;
81d249a78fSBarry Smith 
82d249a78fSBarry Smith   PetscFunctionBegin;
839566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ts->vec_sol, &Y));
849566063dSJacob Faibussowitsch   PetscCall(VecGetSize(ts->vec_sol, &ND));
859566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, ND, NULL, &cvode->work));
869566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, ND, NULL, &cvode->workf));
87d249a78fSBarry Smith 
88d249a78fSBarry Smith   LWORK  = 4 * ND * ND + 12 * ND + 20;
89d249a78fSBarry Smith   LIWORK = 3 * ND + 20;
90d249a78fSBarry Smith 
919566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(LWORK, &WORK, LIWORK, &IWORK));
92d249a78fSBarry Smith 
93d249a78fSBarry Smith   /* C --- PARAMETER IN THE DIFFERENTIAL EQUATION */
94d249a78fSBarry Smith   RPAR = 1.0e-6;
95d249a78fSBarry Smith   /* C --- COMPUTE THE JACOBIAN ANALYTICALLY */
96d249a78fSBarry Smith   IJAC = 1;
97d249a78fSBarry Smith   /* C --- JACOBIAN IS A FULL MATRIX */
98d249a78fSBarry Smith   MLJAC = ND;
99d249a78fSBarry Smith   /* C --- DIFFERENTIAL EQUATION IS IN EXPLICIT FORM*/
100d249a78fSBarry Smith   IMAS = 0;
101d249a78fSBarry Smith   /* C --- OUTPUT ROUTINE IS USED DURING INTEGRATION*/
102d249a78fSBarry Smith   IOUT = 1;
103d249a78fSBarry Smith   /* C --- INITIAL VALUES*/
104d249a78fSBarry Smith   X = ts->ptime;
105d249a78fSBarry Smith   /* C --- ENDPOINT OF INTEGRATION */
106d249a78fSBarry Smith   XEND = ts->max_time;
107d249a78fSBarry Smith   /* C --- REQUIRED TOLERANCE */
1086ffeb343SBarry Smith   RTOL = ts->rtol;
1096ffeb343SBarry Smith   ATOL = ts->atol;
110d249a78fSBarry Smith   ITOL = 0;
111d249a78fSBarry Smith   /* C --- INITIAL STEP SIZE */
11269a506b1SBarry Smith   H = ts->time_step;
113d249a78fSBarry Smith 
114e68ccae8SBarry Smith   /* output MUJAC MLMAS IDID; currently all ignored */
115d249a78fSBarry Smith 
116d249a78fSBarry Smith   radau5_(&ND, FVPOL, &X, Y, &XEND, &H, &RTOL, &ATOL, &ITOL, JVPOL, &IJAC, &MLJAC, &MUJAC, FVPOL, &IMAS, &MLMAS, &MUMAS, SOLOUT, &IOUT, WORK, &LWORK, IWORK, &LIWORK, &RPAR, (void *)ts, &IDID);
117d249a78fSBarry Smith 
1189566063dSJacob Faibussowitsch   PetscCall(PetscFree2(WORK, IWORK));
119d249a78fSBarry Smith   PetscFunctionReturn(0);
120d249a78fSBarry Smith }
121d249a78fSBarry Smith 
122*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDestroy_Radau5(TS ts)
123*d71ae5a4SJacob Faibussowitsch {
124d249a78fSBarry Smith   TS_Radau5 *cvode = (TS_Radau5 *)ts->data;
125d249a78fSBarry Smith 
126d249a78fSBarry Smith   PetscFunctionBegin;
1279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cvode->work));
1289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cvode->workf));
1299566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
130d249a78fSBarry Smith   PetscFunctionReturn(0);
131d249a78fSBarry Smith }
132d249a78fSBarry Smith 
133d249a78fSBarry Smith /*MC
134d249a78fSBarry Smith       TSRADAU5 - ODE solver using the RADAU5 package
135d249a78fSBarry Smith 
13695452b02SPatrick Sanan     Notes:
13795452b02SPatrick Sanan     This uses its own nonlinear solver and dense matrix direct solver so PETSc SNES and KSP options do not apply.
138e68ccae8SBarry Smith            Uses its own time-step adaptivity (but uses the TS rtol and atol, and initial timestep)
139e68ccae8SBarry Smith            Uses its own memory for the dense matrix storage and factorization
140e68ccae8SBarry Smith            Can only handle ODEs of the form \cdot{u} = -F(t,u) + G(t,u)
141d249a78fSBarry Smith 
142d249a78fSBarry Smith     Level: beginner
143d249a78fSBarry Smith 
144db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`
145d249a78fSBarry Smith 
146d249a78fSBarry Smith M*/
147*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Radau5(TS ts)
148*d71ae5a4SJacob Faibussowitsch {
149d249a78fSBarry Smith   TS_Radau5 *cvode;
150d249a78fSBarry Smith 
151d249a78fSBarry Smith   PetscFunctionBegin;
152d249a78fSBarry Smith   ts->ops->destroy       = TSDestroy_Radau5;
153d249a78fSBarry Smith   ts->ops->solve         = TSSolve_Radau5;
154d249a78fSBarry Smith   ts->default_adapt_type = TSADAPTNONE;
155d249a78fSBarry Smith 
1564dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&cvode));
157d249a78fSBarry Smith   ts->data = (void *)cvode;
158d249a78fSBarry Smith   PetscFunctionReturn(0);
159d249a78fSBarry Smith }
160