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 void FVPOL(int *N,double *X,double *Y,double *F,double *RPAR,void *IPAR) 12d249a78fSBarry Smith { 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; 18e68ccae8SBarry Smith PetscErrorCode ierr; 196ffeb343SBarry Smith 20e68ccae8SBarry Smith ierr = VecPlaceArray(cvode->work,Y);CHKERRABORT(PETSC_COMM_SELF,ierr); 21e68ccae8SBarry Smith ierr = VecPlaceArray(cvode->workf,F);CHKERRABORT(PETSC_COMM_SELF,ierr); 226ffeb343SBarry Smith 236ffeb343SBarry Smith /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */ 24e68ccae8SBarry Smith ierr = TSGetDM(ts,&dm);CHKERRABORT(PETSC_COMM_SELF,ierr); 25e68ccae8SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRABORT(PETSC_COMM_SELF,ierr); 26e68ccae8SBarry Smith ierr = DMTSGetIFunction(dm,&ifunction,NULL);CHKERRABORT(PETSC_COMM_SELF,ierr); 276ffeb343SBarry Smith if (!ifunction) { 28e68ccae8SBarry Smith ierr = TSComputeRHSFunction(ts,*X,cvode->work,cvode->workf);CHKERRABORT(PETSC_COMM_SELF,ierr); 29e68ccae8SBarry Smith } else { /* If rhsfunction is also set, this computes both parts and scale them to the right hand side */ 306ffeb343SBarry Smith Vec yydot; 316ffeb343SBarry Smith 32e68ccae8SBarry Smith ierr = VecDuplicate(cvode->work,&yydot);CHKERRABORT(PETSC_COMM_SELF,ierr); 33e68ccae8SBarry Smith ierr = VecZeroEntries(yydot);CHKERRABORT(PETSC_COMM_SELF,ierr); 34e68ccae8SBarry Smith ierr = TSComputeIFunction(ts,*X,cvode->work,yydot,cvode->workf,PETSC_FALSE);CHKERRABORT(PETSC_COMM_SELF,ierr); 35e68ccae8SBarry Smith ierr = VecScale(cvode->workf,-1.);CHKERRABORT(PETSC_COMM_SELF,ierr); 36e68ccae8SBarry Smith ierr = VecDestroy(&yydot);CHKERRABORT(PETSC_COMM_SELF,ierr); 376ffeb343SBarry Smith } 386ffeb343SBarry Smith 39e68ccae8SBarry Smith ierr = VecResetArray(cvode->work);CHKERRABORT(PETSC_COMM_SELF,ierr); 40e68ccae8SBarry Smith ierr = VecResetArray(cvode->workf);CHKERRABORT(PETSC_COMM_SELF,ierr); 416ffeb343SBarry Smith } 426ffeb343SBarry Smith 438acf3572SBarry Smith void JVPOL(PetscInt *N,PetscScalar *X,PetscScalar *Y,PetscScalar *DFY,int *LDFY,PetscScalar *RPAR,void *IPAR) 446ffeb343SBarry Smith { 458acf3572SBarry Smith TS ts = (TS) IPAR; 468acf3572SBarry Smith TS_Radau5 *cvode = (TS_Radau5*)ts->data; 478acf3572SBarry Smith Vec yydot; 488acf3572SBarry Smith Mat mat; 498acf3572SBarry Smith PetscInt n; 50e68ccae8SBarry Smith PetscErrorCode ierr; 51d249a78fSBarry Smith 52e68ccae8SBarry Smith ierr = VecPlaceArray(cvode->work,Y);CHKERRABORT(PETSC_COMM_SELF,ierr); 53e68ccae8SBarry Smith ierr = VecDuplicate(cvode->work,&yydot);CHKERRABORT(PETSC_COMM_SELF,ierr); 54e68ccae8SBarry Smith ierr = VecGetSize(yydot,&n);CHKERRABORT(PETSC_COMM_SELF,ierr); 55e68ccae8SBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,n,n,DFY,&mat);CHKERRABORT(PETSC_COMM_SELF,ierr); 56e68ccae8SBarry Smith ierr = VecZeroEntries(yydot);CHKERRABORT(PETSC_COMM_SELF,ierr); 57e68ccae8SBarry Smith ierr = TSComputeIJacobian(ts,*X,cvode->work,yydot,0,mat,mat,PETSC_FALSE);CHKERRABORT(PETSC_COMM_SELF,ierr); 58e68ccae8SBarry Smith ierr = MatScale(mat,-1.0);CHKERRABORT(PETSC_COMM_SELF,ierr); 59e68ccae8SBarry Smith ierr = MatDestroy(&mat);CHKERRABORT(PETSC_COMM_SELF,ierr); 60e68ccae8SBarry Smith ierr = VecDestroy(&yydot);CHKERRABORT(PETSC_COMM_SELF,ierr); 61e68ccae8SBarry Smith ierr = VecResetArray(cvode->work);CHKERRABORT(PETSC_COMM_SELF,ierr); 628acf3572SBarry Smith } 638acf3572SBarry Smith 64d249a78fSBarry Smith void SOLOUT(int *NR,double *XOLD,double *X, double *Y,double *CONT,double *LRC,int *N,double *RPAR,void *IPAR,int *IRTRN) 65d249a78fSBarry Smith { 66d249a78fSBarry Smith TS ts = (TS) IPAR; 67d249a78fSBarry Smith TS_Radau5 *cvode = (TS_Radau5*)ts->data; 68e68ccae8SBarry Smith PetscErrorCode ierr; 69d249a78fSBarry Smith 70e68ccae8SBarry Smith ierr = VecPlaceArray(cvode->work,Y);CHKERRABORT(PETSC_COMM_SELF,ierr); 71d249a78fSBarry Smith ts->time_step = *X - *XOLD; 72e68ccae8SBarry Smith ierr = TSMonitor(ts,*NR-1,*X,cvode->work);CHKERRABORT(PETSC_COMM_SELF,ierr); 73e68ccae8SBarry Smith ierr = VecResetArray(cvode->work);CHKERRABORT(PETSC_COMM_SELF,ierr); 74d249a78fSBarry Smith } 75d249a78fSBarry Smith 76d249a78fSBarry Smith 77d249a78fSBarry 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*); 78d249a78fSBarry Smith 79d249a78fSBarry Smith PetscErrorCode TSSolve_Radau5(TS ts) 80d249a78fSBarry Smith { 81d249a78fSBarry Smith TS_Radau5 *cvode = (TS_Radau5*)ts->data; 82d249a78fSBarry Smith PetscErrorCode ierr; 83d249a78fSBarry Smith PetscScalar *Y,*WORK,X,XEND,RTOL,ATOL,H,RPAR; 84d249a78fSBarry Smith PetscInt ND,*IWORK,LWORK,LIWORK,MUJAC,MLMAS,MUMAS,IDID,ITOL; 85d249a78fSBarry Smith int IJAC,MLJAC,IMAS,IOUT; 86d249a78fSBarry Smith 87d249a78fSBarry Smith PetscFunctionBegin; 88d249a78fSBarry Smith ierr = VecGetArray(ts->vec_sol,&Y);CHKERRQ(ierr); 89d249a78fSBarry Smith ierr = VecGetSize(ts->vec_sol,&ND);CHKERRQ(ierr); 90d249a78fSBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,ND,NULL,&cvode->work);CHKERRQ(ierr); 916ffeb343SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,ND,NULL,&cvode->workf);CHKERRQ(ierr); 92d249a78fSBarry Smith 93d249a78fSBarry Smith LWORK = 4*ND*ND+12*ND+20; 94d249a78fSBarry Smith LIWORK = 3*ND+20; 95d249a78fSBarry Smith 96d249a78fSBarry Smith ierr = PetscCalloc1(LWORK,&WORK);CHKERRQ(ierr); 97d249a78fSBarry Smith ierr = PetscCalloc1(LIWORK,&IWORK);CHKERRQ(ierr); 98d249a78fSBarry Smith 99d249a78fSBarry Smith /* C --- PARAMETER IN THE DIFFERENTIAL EQUATION */ 100d249a78fSBarry Smith RPAR=1.0e-6; 101d249a78fSBarry Smith /* C --- COMPUTE THE JACOBIAN ANALYTICALLY */ 102d249a78fSBarry Smith IJAC=1; 103d249a78fSBarry Smith /* C --- JACOBIAN IS A FULL MATRIX */ 104d249a78fSBarry Smith MLJAC=ND; 105d249a78fSBarry Smith /* C --- DIFFERENTIAL EQUATION IS IN EXPLICIT FORM*/ 106d249a78fSBarry Smith IMAS=0; 107d249a78fSBarry Smith /* C --- OUTPUT ROUTINE IS USED DURING INTEGRATION*/ 108d249a78fSBarry Smith IOUT=1; 109d249a78fSBarry Smith /* C --- INITIAL VALUES*/ 110d249a78fSBarry Smith X = ts->ptime; 111d249a78fSBarry Smith /* C --- ENDPOINT OF INTEGRATION */ 112d249a78fSBarry Smith XEND = ts->max_time; 113d249a78fSBarry Smith /* C --- REQUIRED TOLERANCE */ 1146ffeb343SBarry Smith RTOL = ts->rtol; 1156ffeb343SBarry Smith ATOL = ts->atol; 116d249a78fSBarry Smith ITOL=0; 117d249a78fSBarry Smith /* C --- INITIAL STEP SIZE */ 11869a506b1SBarry Smith H = ts->time_step; 119d249a78fSBarry Smith 120e68ccae8SBarry Smith /* output MUJAC MLMAS IDID; currently all ignored */ 121d249a78fSBarry Smith 122d249a78fSBarry 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); 123d249a78fSBarry Smith 124d249a78fSBarry Smith ierr = PetscFree(WORK);CHKERRQ(ierr); 125d249a78fSBarry Smith ierr = PetscFree(IWORK);CHKERRQ(ierr); 126d249a78fSBarry Smith PetscFunctionReturn(0); 127d249a78fSBarry Smith } 128d249a78fSBarry Smith 129d249a78fSBarry Smith PetscErrorCode TSDestroy_Radau5(TS ts) 130d249a78fSBarry Smith { 131d249a78fSBarry Smith TS_Radau5 *cvode = (TS_Radau5*)ts->data; 132d249a78fSBarry Smith PetscErrorCode ierr; 133d249a78fSBarry Smith 134d249a78fSBarry Smith PetscFunctionBegin; 135d249a78fSBarry Smith ierr = VecDestroy(&cvode->work);CHKERRQ(ierr); 1366ffeb343SBarry Smith ierr = VecDestroy(&cvode->workf);CHKERRQ(ierr); 137d249a78fSBarry Smith ierr = PetscFree(ts->data);CHKERRQ(ierr); 138d249a78fSBarry Smith PetscFunctionReturn(0); 139d249a78fSBarry Smith } 140d249a78fSBarry Smith 141d249a78fSBarry Smith /*MC 142d249a78fSBarry Smith TSRADAU5 - ODE solver using the RADAU5 package 143d249a78fSBarry Smith 144*95452b02SPatrick Sanan Notes: 145*95452b02SPatrick Sanan This uses its own nonlinear solver and dense matrix direct solver so PETSc SNES and KSP options do not apply. 146e68ccae8SBarry Smith Uses its own time-step adaptivity (but uses the TS rtol and atol, and initial timestep) 147e68ccae8SBarry Smith Uses its own memory for the dense matrix storage and factorization 148e68ccae8SBarry Smith Can only handle ODEs of the form \cdot{u} = -F(t,u) + G(t,u) 149d249a78fSBarry Smith 150d249a78fSBarry Smith Level: beginner 151d249a78fSBarry Smith 152e68ccae8SBarry Smith .seealso: TSCreate(), TS, TSSetType() 153d249a78fSBarry Smith 154d249a78fSBarry Smith M*/ 155d249a78fSBarry Smith PETSC_EXTERN PetscErrorCode TSCreate_Radau5(TS ts) 156d249a78fSBarry Smith { 157d249a78fSBarry Smith TS_Radau5 *cvode; 158d249a78fSBarry Smith PetscErrorCode ierr; 159d249a78fSBarry Smith 160d249a78fSBarry Smith PetscFunctionBegin; 161d249a78fSBarry Smith ts->ops->destroy = TSDestroy_Radau5; 162d249a78fSBarry Smith ts->ops->solve = TSSolve_Radau5; 163d249a78fSBarry Smith ts->default_adapt_type = TSADAPTNONE; 164d249a78fSBarry Smith 165d249a78fSBarry Smith ierr = PetscNewLog(ts,&cvode);CHKERRQ(ierr); 166d249a78fSBarry Smith ts->data = (void*)cvode; 167d249a78fSBarry Smith PetscFunctionReturn(0); 168d249a78fSBarry Smith } 169