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