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 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 { 546ffeb343SBarry Smith TS ts = (TS) IPAR; 556ffeb343SBarry Smith TS_Radau5 *cvode = (TS_Radau5*)ts->data; 566ffeb343SBarry Smith DM dm; 576ffeb343SBarry Smith DMTS tsdm; 586ffeb343SBarry Smith TSIFunction ifunction; 596ffeb343SBarry Smith 606ffeb343SBarry Smith VecPlaceArray(cvode->work,Y); 616ffeb343SBarry Smith VecPlaceArray(cvode->workf,F); 626ffeb343SBarry Smith 636ffeb343SBarry Smith /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */ 646ffeb343SBarry Smith TSGetDM(ts,&dm); 656ffeb343SBarry Smith DMGetDMTS(dm,&tsdm); 666ffeb343SBarry Smith DMTSGetIFunction(dm,&ifunction,NULL); 676ffeb343SBarry Smith if (!ifunction) { 686ffeb343SBarry Smith TSComputeRHSFunction(ts,*X,cvode->work,cvode->workf); 696ffeb343SBarry Smith } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */ 706ffeb343SBarry Smith Vec yydot; 716ffeb343SBarry Smith 726ffeb343SBarry Smith VecDuplicate(cvode->work,&yydot); 736ffeb343SBarry Smith VecZeroEntries(yydot); 746ffeb343SBarry Smith TSComputeIFunction(ts,*X,cvode->work,yydot,cvode->workf,PETSC_FALSE); 756ffeb343SBarry Smith VecScale(cvode->workf,-1.); 766ffeb343SBarry Smith VecDestroy(&yydot); 776ffeb343SBarry Smith } 786ffeb343SBarry Smith 796ffeb343SBarry Smith VecResetArray(cvode->work); 806ffeb343SBarry Smith VecResetArray(cvode->workf); 816ffeb343SBarry Smith } 826ffeb343SBarry Smith 83*8acf3572SBarry Smith void JVPOL(PetscInt *N,PetscScalar *X,PetscScalar *Y,PetscScalar *DFY,int *LDFY,PetscScalar *RPAR,void *IPAR) 846ffeb343SBarry Smith { 85*8acf3572SBarry Smith TS ts = (TS) IPAR; 86*8acf3572SBarry Smith TS_Radau5 *cvode = (TS_Radau5*)ts->data; 87*8acf3572SBarry Smith Vec yydot; 88*8acf3572SBarry Smith Mat mat; 89*8acf3572SBarry Smith PetscInt n; 90d249a78fSBarry Smith 91*8acf3572SBarry Smith VecPlaceArray(cvode->work,Y); 92*8acf3572SBarry Smith VecDuplicate(cvode->work,&yydot); 93*8acf3572SBarry Smith VecGetSize(yydot,&n); 94*8acf3572SBarry Smith MatCreateSeqDense(PETSC_COMM_SELF,n,n,DFY,&mat); 95*8acf3572SBarry Smith VecZeroEntries(yydot); 96*8acf3572SBarry Smith TSComputeIJacobian(ts,*X,cvode->work,yydot,0,mat,mat,PETSC_FALSE); 97*8acf3572SBarry Smith MatScale(mat,-1.0); 98*8acf3572SBarry Smith MatDestroy(&mat); 99*8acf3572SBarry Smith VecDestroy(&yydot); 100*8acf3572SBarry Smith VecResetArray(cvode->work); 101*8acf3572SBarry Smith 102*8acf3572SBarry Smith /* DFY[0]=0.0; 103*8acf3572SBarry Smith DFY[1]=(-2.0*Y[0]*Y[1]-1.0)/(*RPAR); 104*8acf3572SBarry Smith DFY[2]=1.0; 105*8acf3572SBarry Smith DFY[3]=(1.0-Y[0]*Y[0])/(*RPAR); */ 106*8acf3572SBarry Smith } 107*8acf3572SBarry Smith 108*8acf3572SBarry Smith #ifdef foo 109d249a78fSBarry Smith void JVPOL(PetscInt *N,PetscScalar *X,PetscScalar *Y,PetscScalar *DFY,int *LDFY,PetscScalar *RPAR,void *IPAR) 110d249a78fSBarry Smith { 111d249a78fSBarry Smith DFY[0]=0.0; 112d249a78fSBarry Smith DFY[1]=(-2.0*Y[0]*Y[1]-1.0)/(*RPAR); 113d249a78fSBarry Smith DFY[2]=1.0; 114d249a78fSBarry Smith DFY[3]=(1.0-Y[0]*Y[0])/(*RPAR); 115d249a78fSBarry Smith } 116*8acf3572SBarry Smith #endif 117*8acf3572SBarry Smith 118d249a78fSBarry Smith 119d249a78fSBarry Smith void SOLOUT(int *NR,double *XOLD,double *X, double *Y,double *CONT,double *LRC,int *N,double *RPAR,void *IPAR,int *IRTRN) 120d249a78fSBarry Smith { 121d249a78fSBarry Smith TS ts = (TS) IPAR; 122d249a78fSBarry Smith TS_Radau5 *cvode = (TS_Radau5*)ts->data; 123d249a78fSBarry Smith 124d249a78fSBarry Smith VecPlaceArray(cvode->work,Y); 125d249a78fSBarry Smith ts->time_step = *X - *XOLD; 126d249a78fSBarry Smith TSMonitor(ts,*NR-1,*X,cvode->work); 127d249a78fSBarry Smith VecResetArray(cvode->work); 128d249a78fSBarry Smith PetscPrintf(PETSC_COMM_SELF,"X = %g Y = %g %g NSTEP = %d dt = %g\n",*X,Y[0],Y[1],*NR-1,*X - *XOLD); 129d249a78fSBarry Smith } 130d249a78fSBarry Smith 131d249a78fSBarry Smith 132d249a78fSBarry 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*); 133d249a78fSBarry Smith 134d249a78fSBarry Smith PetscErrorCode TSSolve_Radau5(TS ts) 135d249a78fSBarry Smith { 136d249a78fSBarry Smith TS_Radau5 *cvode = (TS_Radau5*)ts->data; 137d249a78fSBarry Smith PetscErrorCode ierr; 138d249a78fSBarry Smith PetscScalar *Y,*WORK,X,XEND,RTOL,ATOL,H,RPAR; 139d249a78fSBarry Smith PetscInt ND,*IWORK,LWORK,LIWORK,MUJAC,MLMAS,MUMAS,IDID,ITOL; 140d249a78fSBarry Smith int IJAC,MLJAC,IMAS,IOUT; 141d249a78fSBarry Smith 142d249a78fSBarry Smith PetscFunctionBegin; 143d249a78fSBarry Smith ierr = VecGetArray(ts->vec_sol,&Y);CHKERRQ(ierr); 144d249a78fSBarry Smith ierr = VecGetSize(ts->vec_sol,&ND);CHKERRQ(ierr); 145d249a78fSBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,ND,NULL,&cvode->work);CHKERRQ(ierr); 1466ffeb343SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,ND,NULL,&cvode->workf);CHKERRQ(ierr); 147d249a78fSBarry Smith 148d249a78fSBarry Smith LWORK = 4*ND*ND+12*ND+20; 149d249a78fSBarry Smith LIWORK = 3*ND+20; 150d249a78fSBarry Smith 151d249a78fSBarry Smith ierr = PetscCalloc1(LWORK,&WORK);CHKERRQ(ierr); 152d249a78fSBarry Smith ierr = PetscCalloc1(LIWORK,&IWORK);CHKERRQ(ierr); 153d249a78fSBarry Smith 154d249a78fSBarry Smith /* C --- PARAMETER IN THE DIFFERENTIAL EQUATION */ 155d249a78fSBarry Smith RPAR=1.0e-6; 156d249a78fSBarry Smith /* C --- COMPUTE THE JACOBIAN ANALYTICALLY */ 157d249a78fSBarry Smith IJAC=1; 158d249a78fSBarry Smith /* C --- JACOBIAN IS A FULL MATRIX */ 159d249a78fSBarry Smith MLJAC=ND; 160d249a78fSBarry Smith /* C --- DIFFERENTIAL EQUATION IS IN EXPLICIT FORM*/ 161d249a78fSBarry Smith IMAS=0; 162d249a78fSBarry Smith /* C --- OUTPUT ROUTINE IS USED DURING INTEGRATION*/ 163d249a78fSBarry Smith IOUT=1; 164d249a78fSBarry Smith /* C --- INITIAL VALUES*/ 165d249a78fSBarry Smith X = ts->ptime; 166d249a78fSBarry Smith /* C --- ENDPOINT OF INTEGRATION */ 167d249a78fSBarry Smith XEND = ts->max_time; 168d249a78fSBarry Smith /* C --- REQUIRED TOLERANCE */ 1696ffeb343SBarry Smith RTOL = ts->rtol; 1706ffeb343SBarry Smith ATOL = ts->atol; 171d249a78fSBarry Smith ITOL=0; 172d249a78fSBarry Smith /* C --- INITIAL STEP SIZE */ 17369a506b1SBarry Smith H = ts->time_step; 174d249a78fSBarry Smith 175d249a78fSBarry Smith /* output MUJAC MLMAS IDID */ 176d249a78fSBarry Smith 177d249a78fSBarry Smith CHKMEMQ; 178d249a78fSBarry 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); 179d249a78fSBarry Smith CHKMEMQ; 180d249a78fSBarry Smith 181d249a78fSBarry Smith ierr = PetscFree(WORK);CHKERRQ(ierr); 182d249a78fSBarry Smith ierr = PetscFree(IWORK);CHKERRQ(ierr); 183d249a78fSBarry Smith PetscFunctionReturn(0); 184d249a78fSBarry Smith } 185d249a78fSBarry Smith 186d249a78fSBarry Smith PetscErrorCode TSDestroy_Radau5(TS ts) 187d249a78fSBarry Smith { 188d249a78fSBarry Smith TS_Radau5 *cvode = (TS_Radau5*)ts->data; 189d249a78fSBarry Smith PetscErrorCode ierr; 190d249a78fSBarry Smith 191d249a78fSBarry Smith PetscFunctionBegin; 192d249a78fSBarry Smith ierr = VecDestroy(&cvode->work);CHKERRQ(ierr); 1936ffeb343SBarry Smith ierr = VecDestroy(&cvode->workf);CHKERRQ(ierr); 194d249a78fSBarry Smith ierr = PetscFree(ts->data);CHKERRQ(ierr); 195d249a78fSBarry Smith PetscFunctionReturn(0); 196d249a78fSBarry Smith } 197d249a78fSBarry Smith 198d249a78fSBarry Smith PetscErrorCode TSSetFromOptions_Radau5(PetscOptionItems *PetscOptionsObject,TS ts) 199d249a78fSBarry Smith { 200d249a78fSBarry Smith TS_Radau5 *cvode = (TS_Radau5*)ts->data; 201d249a78fSBarry Smith PetscErrorCode ierr; 202d249a78fSBarry Smith 203d249a78fSBarry Smith PetscFunctionBegin; 204d249a78fSBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"RADAU5 ODE solver options");CHKERRQ(ierr); 205d249a78fSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 206d249a78fSBarry Smith PetscFunctionReturn(0); 207d249a78fSBarry Smith } 208d249a78fSBarry Smith 209d249a78fSBarry Smith PetscErrorCode TSView_Radau5(TS ts,PetscViewer viewer) 210d249a78fSBarry Smith { 211d249a78fSBarry Smith TS_Radau5 *cvode = (TS_Radau5*)ts->data; 212d249a78fSBarry Smith PetscErrorCode ierr; 213d249a78fSBarry Smith 214d249a78fSBarry Smith PetscFunctionBegin; 215d249a78fSBarry Smith PetscFunctionReturn(0); 216d249a78fSBarry Smith } 217d249a78fSBarry Smith 218d249a78fSBarry Smith /*MC 219d249a78fSBarry Smith TSRADAU5 - ODE solver using the RADAU5 package 220d249a78fSBarry Smith 221d249a78fSBarry Smith Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply 222d249a78fSBarry Smith 223d249a78fSBarry Smith Level: beginner 224d249a78fSBarry Smith 225d249a78fSBarry Smith .seealso: TSCreate(), TS, TSSetType(), TSSetExactFinalTime() 226d249a78fSBarry Smith 227d249a78fSBarry Smith M*/ 228d249a78fSBarry Smith PETSC_EXTERN PetscErrorCode TSCreate_Radau5(TS ts) 229d249a78fSBarry Smith { 230d249a78fSBarry Smith TS_Radau5 *cvode; 231d249a78fSBarry Smith PetscErrorCode ierr; 232d249a78fSBarry Smith 233d249a78fSBarry Smith PetscFunctionBegin; 234d249a78fSBarry Smith ts->ops->destroy = TSDestroy_Radau5; 235d249a78fSBarry Smith ts->ops->view = TSView_Radau5; 236d249a78fSBarry Smith ts->ops->solve = TSSolve_Radau5; 237d249a78fSBarry Smith ts->ops->setfromoptions = TSSetFromOptions_Radau5; 238d249a78fSBarry Smith ts->default_adapt_type = TSADAPTNONE; 239d249a78fSBarry Smith 240d249a78fSBarry Smith ierr = PetscNewLog(ts,&cvode);CHKERRQ(ierr); 241d249a78fSBarry Smith ts->data = (void*)cvode; 242d249a78fSBarry Smith PetscFunctionReturn(0); 243d249a78fSBarry Smith } 244