xref: /petsc/src/ts/impls/implicit/radau5/radau5.c (revision 6ffeb3436a504f614c4b292f2b28b43da6b401b2)
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 {
8*6ffeb343SBarry Smith   Vec work,workf;
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 {
54*6ffeb343SBarry Smith   TS          ts = (TS) IPAR;
55*6ffeb343SBarry Smith   TS_Radau5   *cvode = (TS_Radau5*)ts->data;
56*6ffeb343SBarry Smith   DM          dm;
57*6ffeb343SBarry Smith   DMTS        tsdm;
58*6ffeb343SBarry Smith   TSIFunction ifunction;
59*6ffeb343SBarry Smith 
60*6ffeb343SBarry Smith   VecPlaceArray(cvode->work,Y);
61*6ffeb343SBarry Smith   VecPlaceArray(cvode->workf,F);
62*6ffeb343SBarry Smith 
63*6ffeb343SBarry Smith   /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */
64*6ffeb343SBarry Smith   TSGetDM(ts,&dm);
65*6ffeb343SBarry Smith   DMGetDMTS(dm,&tsdm);
66*6ffeb343SBarry Smith   DMTSGetIFunction(dm,&ifunction,NULL);
67*6ffeb343SBarry Smith   if (!ifunction) {
68*6ffeb343SBarry Smith     TSComputeRHSFunction(ts,*X,cvode->work,cvode->workf);
69*6ffeb343SBarry Smith   } else {                      /* If rhsfunction is also set, this computes both parts and shifts them to the right */
70*6ffeb343SBarry Smith     Vec yydot;
71*6ffeb343SBarry Smith 
72*6ffeb343SBarry Smith     VecDuplicate(cvode->work,&yydot);
73*6ffeb343SBarry Smith     VecZeroEntries(yydot);
74*6ffeb343SBarry Smith     TSComputeIFunction(ts,*X,cvode->work,yydot,cvode->workf,PETSC_FALSE);
75*6ffeb343SBarry Smith     VecScale(cvode->workf,-1.);
76*6ffeb343SBarry Smith     VecDestroy(&yydot);
77*6ffeb343SBarry Smith   }
78*6ffeb343SBarry Smith 
79*6ffeb343SBarry Smith   VecResetArray(cvode->work);
80*6ffeb343SBarry Smith   VecResetArray(cvode->workf);
81*6ffeb343SBarry Smith }
82*6ffeb343SBarry Smith 
83*6ffeb343SBarry Smith #ifdef foo
84*6ffeb343SBarry Smith void FVPOL(int *N,double *X,double *Y,double *F,double *RPAR,void *IPAR)
85*6ffeb343SBarry Smith {
86d249a78fSBarry Smith   F[0]=Y[1];
87d249a78fSBarry Smith   F[1]=((1-Y[0]*Y[0])*Y[1]-Y[0])/(*RPAR);
88d249a78fSBarry Smith }
89*6ffeb343SBarry Smith #endif
90d249a78fSBarry Smith 
91d249a78fSBarry Smith void JVPOL(PetscInt *N,PetscScalar *X,PetscScalar *Y,PetscScalar *DFY,int *LDFY,PetscScalar *RPAR,void *IPAR)
92d249a78fSBarry Smith {
93d249a78fSBarry Smith   DFY[0]=0.0;
94d249a78fSBarry Smith   DFY[1]=(-2.0*Y[0]*Y[1]-1.0)/(*RPAR);
95d249a78fSBarry Smith   DFY[2]=1.0;
96d249a78fSBarry Smith   DFY[3]=(1.0-Y[0]*Y[0])/(*RPAR);
97d249a78fSBarry Smith }
98d249a78fSBarry Smith 
99d249a78fSBarry Smith void SOLOUT(int *NR,double *XOLD,double *X, double *Y,double *CONT,double *LRC,int *N,double *RPAR,void *IPAR,int *IRTRN)
100d249a78fSBarry Smith {
101d249a78fSBarry Smith   TS        ts = (TS) IPAR;
102d249a78fSBarry Smith   TS_Radau5 *cvode = (TS_Radau5*)ts->data;
103d249a78fSBarry Smith 
104d249a78fSBarry Smith   VecPlaceArray(cvode->work,Y);
105d249a78fSBarry Smith   ts->time_step = *X - *XOLD;
106d249a78fSBarry Smith   TSMonitor(ts,*NR-1,*X,cvode->work);
107d249a78fSBarry Smith   VecResetArray(cvode->work);
108d249a78fSBarry Smith   PetscPrintf(PETSC_COMM_SELF,"X = %g Y = %g %g NSTEP = %d dt = %g\n",*X,Y[0],Y[1],*NR-1,*X - *XOLD);
109d249a78fSBarry Smith }
110d249a78fSBarry Smith 
111d249a78fSBarry Smith 
112d249a78fSBarry 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*);
113d249a78fSBarry Smith 
114d249a78fSBarry Smith PetscErrorCode TSSolve_Radau5(TS ts)
115d249a78fSBarry Smith {
116d249a78fSBarry Smith   TS_Radau5      *cvode = (TS_Radau5*)ts->data;
117d249a78fSBarry Smith   PetscErrorCode ierr;
118d249a78fSBarry Smith   PetscScalar    *Y,*WORK,X,XEND,RTOL,ATOL,H,RPAR;
119d249a78fSBarry Smith   PetscInt       ND,*IWORK,LWORK,LIWORK,MUJAC,MLMAS,MUMAS,IDID,ITOL;
120d249a78fSBarry Smith   int            IJAC,MLJAC,IMAS,IOUT;
121d249a78fSBarry Smith 
122d249a78fSBarry Smith   PetscFunctionBegin;
123d249a78fSBarry Smith   ierr = VecGetArray(ts->vec_sol,&Y);CHKERRQ(ierr);
124d249a78fSBarry Smith   ierr = VecGetSize(ts->vec_sol,&ND);CHKERRQ(ierr);
125d249a78fSBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,ND,NULL,&cvode->work);CHKERRQ(ierr);
126*6ffeb343SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,ND,NULL,&cvode->workf);CHKERRQ(ierr);
127d249a78fSBarry Smith 
128d249a78fSBarry Smith   LWORK  = 4*ND*ND+12*ND+20;
129d249a78fSBarry Smith   LIWORK = 3*ND+20;
130d249a78fSBarry Smith 
131d249a78fSBarry Smith   ierr = PetscCalloc1(LWORK,&WORK);CHKERRQ(ierr);
132d249a78fSBarry Smith   ierr = PetscCalloc1(LIWORK,&IWORK);CHKERRQ(ierr);
133d249a78fSBarry Smith 
134d249a78fSBarry Smith   /* C --- PARAMETER IN THE DIFFERENTIAL EQUATION */
135d249a78fSBarry Smith   RPAR=1.0e-6;
136d249a78fSBarry Smith   /* C --- COMPUTE THE JACOBIAN ANALYTICALLY */
137d249a78fSBarry Smith   IJAC=1;
138d249a78fSBarry Smith   /* C --- JACOBIAN IS A FULL MATRIX */
139d249a78fSBarry Smith   MLJAC=ND;
140d249a78fSBarry Smith   /* C --- DIFFERENTIAL EQUATION IS IN EXPLICIT FORM*/
141d249a78fSBarry Smith   IMAS=0;
142d249a78fSBarry Smith   /* C --- OUTPUT ROUTINE IS USED DURING INTEGRATION*/
143d249a78fSBarry Smith   IOUT=1;
144d249a78fSBarry Smith   /* C --- INITIAL VALUES*/
145d249a78fSBarry Smith   X = ts->ptime;
146d249a78fSBarry Smith   /* C --- ENDPOINT OF INTEGRATION */
147d249a78fSBarry Smith   XEND = ts->max_time;
148d249a78fSBarry Smith   /* C --- REQUIRED TOLERANCE */
149*6ffeb343SBarry Smith   RTOL = ts->rtol;
150*6ffeb343SBarry Smith   ATOL = ts->atol;
151d249a78fSBarry Smith   ITOL=0;
152d249a78fSBarry Smith   /* C --- INITIAL STEP SIZE */
15369a506b1SBarry Smith   H = ts->time_step;
154d249a78fSBarry Smith 
155d249a78fSBarry Smith   /* output MUJAC MLMAS IDID */
156d249a78fSBarry Smith 
157d249a78fSBarry Smith   CHKMEMQ;
158d249a78fSBarry 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);
159d249a78fSBarry Smith   CHKMEMQ;
160d249a78fSBarry Smith 
161d249a78fSBarry Smith   ierr = PetscFree(WORK);CHKERRQ(ierr);
162d249a78fSBarry Smith   ierr = PetscFree(IWORK);CHKERRQ(ierr);
163d249a78fSBarry Smith   PetscFunctionReturn(0);
164d249a78fSBarry Smith }
165d249a78fSBarry Smith 
166d249a78fSBarry Smith PetscErrorCode TSDestroy_Radau5(TS ts)
167d249a78fSBarry Smith {
168d249a78fSBarry Smith   TS_Radau5      *cvode = (TS_Radau5*)ts->data;
169d249a78fSBarry Smith   PetscErrorCode ierr;
170d249a78fSBarry Smith 
171d249a78fSBarry Smith   PetscFunctionBegin;
172d249a78fSBarry Smith   ierr = VecDestroy(&cvode->work);CHKERRQ(ierr);
173*6ffeb343SBarry Smith   ierr = VecDestroy(&cvode->workf);CHKERRQ(ierr);
174d249a78fSBarry Smith   ierr = PetscFree(ts->data);CHKERRQ(ierr);
175d249a78fSBarry Smith   PetscFunctionReturn(0);
176d249a78fSBarry Smith }
177d249a78fSBarry Smith 
178d249a78fSBarry Smith PetscErrorCode TSSetFromOptions_Radau5(PetscOptionItems *PetscOptionsObject,TS ts)
179d249a78fSBarry Smith {
180d249a78fSBarry Smith   TS_Radau5      *cvode = (TS_Radau5*)ts->data;
181d249a78fSBarry Smith   PetscErrorCode ierr;
182d249a78fSBarry Smith 
183d249a78fSBarry Smith   PetscFunctionBegin;
184d249a78fSBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RADAU5 ODE solver options");CHKERRQ(ierr);
185d249a78fSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
186d249a78fSBarry Smith   PetscFunctionReturn(0);
187d249a78fSBarry Smith }
188d249a78fSBarry Smith 
189d249a78fSBarry Smith PetscErrorCode TSView_Radau5(TS ts,PetscViewer viewer)
190d249a78fSBarry Smith {
191d249a78fSBarry Smith   TS_Radau5      *cvode = (TS_Radau5*)ts->data;
192d249a78fSBarry Smith   PetscErrorCode ierr;
193d249a78fSBarry Smith 
194d249a78fSBarry Smith   PetscFunctionBegin;
195d249a78fSBarry Smith   PetscFunctionReturn(0);
196d249a78fSBarry Smith }
197d249a78fSBarry Smith 
198d249a78fSBarry Smith /*MC
199d249a78fSBarry Smith       TSRADAU5 - ODE solver using the RADAU5 package
200d249a78fSBarry Smith 
201d249a78fSBarry Smith     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
202d249a78fSBarry Smith 
203d249a78fSBarry Smith     Level: beginner
204d249a78fSBarry Smith 
205d249a78fSBarry Smith .seealso:  TSCreate(), TS, TSSetType(), TSSetExactFinalTime()
206d249a78fSBarry Smith 
207d249a78fSBarry Smith M*/
208d249a78fSBarry Smith PETSC_EXTERN PetscErrorCode TSCreate_Radau5(TS ts)
209d249a78fSBarry Smith {
210d249a78fSBarry Smith   TS_Radau5      *cvode;
211d249a78fSBarry Smith   PetscErrorCode ierr;
212d249a78fSBarry Smith 
213d249a78fSBarry Smith   PetscFunctionBegin;
214d249a78fSBarry Smith   ts->ops->destroy        = TSDestroy_Radau5;
215d249a78fSBarry Smith   ts->ops->view           = TSView_Radau5;
216d249a78fSBarry Smith   ts->ops->solve          = TSSolve_Radau5;
217d249a78fSBarry Smith   ts->ops->setfromoptions = TSSetFromOptions_Radau5;
218d249a78fSBarry Smith   ts->default_adapt_type  = TSADAPTNONE;
219d249a78fSBarry Smith 
220d249a78fSBarry Smith   ierr = PetscNewLog(ts,&cvode);CHKERRQ(ierr);
221d249a78fSBarry Smith   ts->data = (void*)cvode;
222d249a78fSBarry Smith   PetscFunctionReturn(0);
223d249a78fSBarry Smith }
224