xref: /petsc/src/ts/impls/implicit/radau5/radau5.c (revision 69a506b14013087b93726df018e15ded3faff328)
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 {
8d249a78fSBarry Smith   Vec work;
9d249a78fSBarry Smith } TS_Radau5;
10d249a78fSBarry Smith 
11d249a78fSBarry Smith #ifdef foo
12d249a78fSBarry Smith /*
13d249a78fSBarry Smith         TSFunction_Radau5 - routine that we provide to Radau5 that applies the right hand side.
14d249a78fSBarry Smith */
15d249a78fSBarry Smith int TSFunction_Radau5(realtype t,N_Vector y,N_Vector ydot,void *ctx)
16d249a78fSBarry Smith {
17d249a78fSBarry Smith   TS             ts = (TS) ctx;
18d249a78fSBarry Smith   DM             dm;
19d249a78fSBarry Smith   DMTS           tsdm;
20d249a78fSBarry Smith   TSIFunction    ifunction;
21d249a78fSBarry Smith   MPI_Comm       comm;
22d249a78fSBarry Smith   TS_Radau5    *cvode = (TS_Radau5*)ts->data;
23d249a78fSBarry Smith   Vec            yy     = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot;
24d249a78fSBarry Smith   PetscScalar    *y_data,*ydot_data;
25d249a78fSBarry Smith   PetscErrorCode ierr;
26d249a78fSBarry Smith 
27d249a78fSBarry Smith   PetscFunctionBegin;
28d249a78fSBarry Smith   ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
29d249a78fSBarry Smith   /* Make the PETSc work vectors yy and yyd point to the arrays in the RADAU5 vectors y and ydot respectively*/
30d249a78fSBarry Smith   y_data    = (PetscScalar*) N_VGetArrayPointer(y);
31d249a78fSBarry Smith   ydot_data = (PetscScalar*) N_VGetArrayPointer(ydot);
32d249a78fSBarry Smith   ierr      = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
33d249a78fSBarry Smith   ierr      = VecPlaceArray(yyd,ydot_data);CHKERRABORT(comm,ierr);
34d249a78fSBarry Smith 
35d249a78fSBarry Smith   /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */
36d249a78fSBarry Smith   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
37d249a78fSBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
38d249a78fSBarry Smith   ierr = DMTSGetIFunction(dm,&ifunction,NULL);CHKERRQ(ierr);
39d249a78fSBarry Smith   if (!ifunction) {
40d249a78fSBarry Smith     ierr = TSComputeRHSFunction(ts,t,yy,yyd);CHKERRQ(ierr);
41d249a78fSBarry Smith   } else {                      /* If rhsfunction is also set, this computes both parts and shifts them to the right */
42d249a78fSBarry Smith     ierr = VecZeroEntries(yydot);CHKERRQ(ierr);
43d249a78fSBarry Smith     ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE);CHKERRABORT(comm,ierr);
44d249a78fSBarry Smith     ierr = VecScale(yyd,-1.);CHKERRQ(ierr);
45d249a78fSBarry Smith   }
46d249a78fSBarry Smith   ierr = VecResetArray(yy);CHKERRABORT(comm,ierr);
47d249a78fSBarry Smith   ierr = VecResetArray(yyd);CHKERRABORT(comm,ierr);
48d249a78fSBarry Smith   PetscFunctionReturn(0);
49d249a78fSBarry Smith }
50d249a78fSBarry Smith #endif
51d249a78fSBarry Smith 
52d249a78fSBarry Smith void FVPOL(int *N,double *X,double *Y,double *F,double *RPAR,void *IPAR)
53d249a78fSBarry Smith {
54d249a78fSBarry Smith   F[0]=Y[1];
55d249a78fSBarry Smith   F[1]=((1-Y[0]*Y[0])*Y[1]-Y[0])/(*RPAR);
56d249a78fSBarry Smith }
57d249a78fSBarry Smith 
58d249a78fSBarry Smith void JVPOL(PetscInt *N,PetscScalar *X,PetscScalar *Y,PetscScalar *DFY,int *LDFY,PetscScalar *RPAR,void *IPAR)
59d249a78fSBarry Smith {
60d249a78fSBarry Smith   DFY[0]=0.0;
61d249a78fSBarry Smith   DFY[1]=(-2.0*Y[0]*Y[1]-1.0)/(*RPAR);
62d249a78fSBarry Smith   DFY[2]=1.0;
63d249a78fSBarry Smith   DFY[3]=(1.0-Y[0]*Y[0])/(*RPAR);
64d249a78fSBarry Smith }
65d249a78fSBarry Smith 
66d249a78fSBarry Smith void SOLOUT(int *NR,double *XOLD,double *X, double *Y,double *CONT,double *LRC,int *N,double *RPAR,void *IPAR,int *IRTRN)
67d249a78fSBarry Smith {
68d249a78fSBarry Smith   TS        ts = (TS) IPAR;
69d249a78fSBarry Smith   TS_Radau5 *cvode = (TS_Radau5*)ts->data;
70d249a78fSBarry Smith 
71d249a78fSBarry Smith   VecPlaceArray(cvode->work,Y);
72d249a78fSBarry Smith   ts->time_step = *X - *XOLD;
73d249a78fSBarry Smith   TSMonitor(ts,*NR-1,*X,cvode->work);
74d249a78fSBarry Smith   VecResetArray(cvode->work);
75d249a78fSBarry Smith   PetscPrintf(PETSC_COMM_SELF,"X = %g Y = %g %g NSTEP = %d dt = %g\n",*X,Y[0],Y[1],*NR-1,*X - *XOLD);
76d249a78fSBarry Smith }
77d249a78fSBarry Smith 
78d249a78fSBarry Smith 
79d249a78fSBarry 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*);
80d249a78fSBarry Smith 
81d249a78fSBarry Smith PetscErrorCode TSSolve_Radau5(TS ts)
82d249a78fSBarry Smith {
83d249a78fSBarry Smith   TS_Radau5      *cvode = (TS_Radau5*)ts->data;
84d249a78fSBarry Smith   PetscErrorCode ierr;
85d249a78fSBarry Smith   PetscScalar    *Y,*WORK,X,XEND,RTOL,ATOL,H,RPAR;
86d249a78fSBarry Smith   PetscInt       ND,*IWORK,LWORK,LIWORK,MUJAC,MLMAS,MUMAS,IDID,ITOL;
87d249a78fSBarry Smith   int            IJAC,MLJAC,IMAS,IOUT;
88d249a78fSBarry Smith 
89d249a78fSBarry Smith   PetscFunctionBegin;
90d249a78fSBarry Smith   ierr = VecGetArray(ts->vec_sol,&Y);CHKERRQ(ierr);
91d249a78fSBarry Smith   ierr = VecGetSize(ts->vec_sol,&ND);CHKERRQ(ierr);
92d249a78fSBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,ND,NULL,&cvode->work);CHKERRQ(ierr);
93d249a78fSBarry Smith 
94d249a78fSBarry Smith   LWORK  = 4*ND*ND+12*ND+20;
95d249a78fSBarry Smith   LIWORK = 3*ND+20;
96d249a78fSBarry Smith 
97d249a78fSBarry Smith   ierr = PetscCalloc1(LWORK,&WORK);CHKERRQ(ierr);
98d249a78fSBarry Smith   ierr = PetscCalloc1(LIWORK,&IWORK);CHKERRQ(ierr);
99d249a78fSBarry Smith 
100d249a78fSBarry Smith   /* C --- PARAMETER IN THE DIFFERENTIAL EQUATION */
101d249a78fSBarry Smith   RPAR=1.0e-6;
102d249a78fSBarry Smith   /* C --- COMPUTE THE JACOBIAN ANALYTICALLY */
103d249a78fSBarry Smith   IJAC=1;
104d249a78fSBarry Smith   /* C --- JACOBIAN IS A FULL MATRIX */
105d249a78fSBarry Smith   MLJAC=ND;
106d249a78fSBarry Smith   /* C --- DIFFERENTIAL EQUATION IS IN EXPLICIT FORM*/
107d249a78fSBarry Smith   IMAS=0;
108d249a78fSBarry Smith   /* C --- OUTPUT ROUTINE IS USED DURING INTEGRATION*/
109d249a78fSBarry Smith   IOUT=1;
110d249a78fSBarry Smith   /* C --- INITIAL VALUES*/
111d249a78fSBarry Smith   X = ts->ptime;
112d249a78fSBarry Smith   /* C --- ENDPOINT OF INTEGRATION */
113d249a78fSBarry Smith   XEND = ts->max_time;
114d249a78fSBarry Smith   /* C --- REQUIRED TOLERANCE */
115d249a78fSBarry Smith   RTOL=1.0e-4;
116d249a78fSBarry Smith   ATOL=1.0*RTOL;
117d249a78fSBarry Smith   ITOL=0;
118d249a78fSBarry Smith   /* C --- INITIAL STEP SIZE */
119*69a506b1SBarry Smith   H = ts->time_step;
120d249a78fSBarry Smith 
121d249a78fSBarry Smith   /* output MUJAC MLMAS IDID */
122d249a78fSBarry Smith 
123d249a78fSBarry Smith   CHKMEMQ;
124d249a78fSBarry 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);
125d249a78fSBarry Smith   CHKMEMQ;
126d249a78fSBarry Smith 
127d249a78fSBarry Smith   ierr = PetscFree(WORK);CHKERRQ(ierr);
128d249a78fSBarry Smith   ierr = PetscFree(IWORK);CHKERRQ(ierr);
129d249a78fSBarry Smith   PetscFunctionReturn(0);
130d249a78fSBarry Smith }
131d249a78fSBarry Smith 
132d249a78fSBarry Smith PetscErrorCode TSDestroy_Radau5(TS ts)
133d249a78fSBarry Smith {
134d249a78fSBarry Smith   TS_Radau5      *cvode = (TS_Radau5*)ts->data;
135d249a78fSBarry Smith   PetscErrorCode ierr;
136d249a78fSBarry Smith 
137d249a78fSBarry Smith   PetscFunctionBegin;
138d249a78fSBarry Smith   ierr = VecDestroy(&cvode->work);CHKERRQ(ierr);
139d249a78fSBarry Smith   ierr = PetscFree(ts->data);CHKERRQ(ierr);
140d249a78fSBarry Smith   PetscFunctionReturn(0);
141d249a78fSBarry Smith }
142d249a78fSBarry Smith 
143d249a78fSBarry Smith PetscErrorCode TSSetFromOptions_Radau5(PetscOptionItems *PetscOptionsObject,TS ts)
144d249a78fSBarry Smith {
145d249a78fSBarry Smith   TS_Radau5      *cvode = (TS_Radau5*)ts->data;
146d249a78fSBarry Smith   PetscErrorCode ierr;
147d249a78fSBarry Smith 
148d249a78fSBarry Smith   PetscFunctionBegin;
149d249a78fSBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RADAU5 ODE solver options");CHKERRQ(ierr);
150d249a78fSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
151d249a78fSBarry Smith   PetscFunctionReturn(0);
152d249a78fSBarry Smith }
153d249a78fSBarry Smith 
154d249a78fSBarry Smith PetscErrorCode TSView_Radau5(TS ts,PetscViewer viewer)
155d249a78fSBarry Smith {
156d249a78fSBarry Smith   TS_Radau5      *cvode = (TS_Radau5*)ts->data;
157d249a78fSBarry Smith   PetscErrorCode ierr;
158d249a78fSBarry Smith 
159d249a78fSBarry Smith   PetscFunctionBegin;
160d249a78fSBarry Smith   PetscFunctionReturn(0);
161d249a78fSBarry Smith }
162d249a78fSBarry Smith 
163d249a78fSBarry Smith /*MC
164d249a78fSBarry Smith       TSRADAU5 - ODE solver using the RADAU5 package
165d249a78fSBarry Smith 
166d249a78fSBarry Smith     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
167d249a78fSBarry Smith 
168d249a78fSBarry Smith     Level: beginner
169d249a78fSBarry Smith 
170d249a78fSBarry Smith .seealso:  TSCreate(), TS, TSSetType(), TSSetExactFinalTime()
171d249a78fSBarry Smith 
172d249a78fSBarry Smith M*/
173d249a78fSBarry Smith PETSC_EXTERN PetscErrorCode TSCreate_Radau5(TS ts)
174d249a78fSBarry Smith {
175d249a78fSBarry Smith   TS_Radau5      *cvode;
176d249a78fSBarry Smith   PetscErrorCode ierr;
177d249a78fSBarry Smith 
178d249a78fSBarry Smith   PetscFunctionBegin;
179d249a78fSBarry Smith   ts->ops->destroy        = TSDestroy_Radau5;
180d249a78fSBarry Smith   ts->ops->view           = TSView_Radau5;
181d249a78fSBarry Smith   ts->ops->solve          = TSSolve_Radau5;
182d249a78fSBarry Smith   ts->ops->setfromoptions = TSSetFromOptions_Radau5;
183d249a78fSBarry Smith   ts->default_adapt_type  = TSADAPTNONE;
184d249a78fSBarry Smith 
185d249a78fSBarry Smith   ierr = PetscNewLog(ts,&cvode);CHKERRQ(ierr);
186d249a78fSBarry Smith   ts->data = (void*)cvode;
187d249a78fSBarry Smith   PetscFunctionReturn(0);
188d249a78fSBarry Smith }
189