17cb94ee6SHong Zhang #define PETSCTS_DLL 27cb94ee6SHong Zhang 37cb94ee6SHong Zhang /* 4db268bfaSHong Zhang Provides a PETSc interface to SUNDIALS/CVODE solver. 57cb94ee6SHong Zhang The interface to PVODE (old version of CVODE) was originally contributed 67cb94ee6SHong Zhang by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik. 728adb3f7SHong Zhang 828adb3f7SHong Zhang Reference: sundials-2.4.0/examples/cvode/parallel/cvDiurnal_kry_p.c 97cb94ee6SHong Zhang */ 1056a740aaSHong Zhang #include "sundials.h" /*I "petscts.h" I*/ 117cb94ee6SHong Zhang 127cb94ee6SHong Zhang /* 137cb94ee6SHong Zhang TSPrecond_Sundials - function that we provide to SUNDIALS to 147cb94ee6SHong Zhang evaluate the preconditioner. 157cb94ee6SHong Zhang */ 167cb94ee6SHong Zhang #undef __FUNCT__ 177cb94ee6SHong Zhang #define __FUNCT__ "TSPrecond_Sundials" 187cb94ee6SHong Zhang PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy, 197cb94ee6SHong Zhang booleantype jok,booleantype *jcurPtr, 207cb94ee6SHong Zhang realtype _gamma,void *P_data, 217cb94ee6SHong Zhang N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3) 227cb94ee6SHong Zhang { 237cb94ee6SHong Zhang TS ts = (TS) P_data; 247cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 257cb94ee6SHong Zhang PC pc = cvode->pc; 267cb94ee6SHong Zhang PetscErrorCode ierr; 277cb94ee6SHong Zhang Mat Jac = ts->B; 287cb94ee6SHong Zhang Vec yy = cvode->w1; 297cb94ee6SHong Zhang PetscScalar one = 1.0,gm; 307cb94ee6SHong Zhang MatStructure str = DIFFERENT_NONZERO_PATTERN; 317cb94ee6SHong Zhang PetscScalar *y_data; 327cb94ee6SHong Zhang 337cb94ee6SHong Zhang PetscFunctionBegin; 347cb94ee6SHong Zhang /* This allows us to construct preconditioners in-place if we like */ 357cb94ee6SHong Zhang ierr = MatSetUnfactored(Jac);CHKERRQ(ierr); 367cb94ee6SHong Zhang 377cb94ee6SHong Zhang /* jok - TRUE means reuse current Jacobian else recompute Jacobian */ 387cb94ee6SHong Zhang if (jok) { 397cb94ee6SHong Zhang ierr = MatCopy(cvode->pmat,Jac,str);CHKERRQ(ierr); 407cb94ee6SHong Zhang *jcurPtr = FALSE; 417cb94ee6SHong Zhang } else { 427cb94ee6SHong Zhang /* make PETSc vector yy point to SUNDIALS vector y */ 437cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(y); 447cb94ee6SHong Zhang ierr = VecPlaceArray(yy,y_data); CHKERRQ(ierr); 454ac7836bSHong Zhang 467cb94ee6SHong Zhang /* compute the Jacobian */ 477cb94ee6SHong Zhang ierr = TSComputeRHSJacobian(ts,ts->ptime,yy,&Jac,&Jac,&str);CHKERRQ(ierr); 487cb94ee6SHong Zhang ierr = VecResetArray(yy); CHKERRQ(ierr); 494ac7836bSHong Zhang 507cb94ee6SHong Zhang /* copy the Jacobian matrix */ 517cb94ee6SHong Zhang if (!cvode->pmat) { 527cb94ee6SHong Zhang ierr = MatDuplicate(Jac,MAT_COPY_VALUES,&cvode->pmat);CHKERRQ(ierr); 537cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->pmat);CHKERRQ(ierr); 542c823083SHong Zhang } else { 557cb94ee6SHong Zhang ierr = MatCopy(Jac,cvode->pmat,str);CHKERRQ(ierr); 567cb94ee6SHong Zhang } 577cb94ee6SHong Zhang *jcurPtr = TRUE; 587cb94ee6SHong Zhang } 597cb94ee6SHong Zhang 607cb94ee6SHong Zhang /* construct I-gamma*Jac */ 617cb94ee6SHong Zhang gm = -_gamma; 627cb94ee6SHong Zhang ierr = MatScale(Jac,gm);CHKERRQ(ierr); 637cb94ee6SHong Zhang ierr = MatShift(Jac,one);CHKERRQ(ierr); 647cb94ee6SHong Zhang 657cb94ee6SHong Zhang ierr = PCSetOperators(pc,Jac,Jac,str);CHKERRQ(ierr); 667cb94ee6SHong Zhang PetscFunctionReturn(0); 677cb94ee6SHong Zhang } 687cb94ee6SHong Zhang 697cb94ee6SHong Zhang /* 707cb94ee6SHong Zhang TSPSolve_Sundials - routine that we provide to Sundials that applies the preconditioner. 717cb94ee6SHong Zhang */ 727cb94ee6SHong Zhang #undef __FUNCT__ 737cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials" 744ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z, 754ac7836bSHong Zhang realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp) 767cb94ee6SHong Zhang { 777cb94ee6SHong Zhang TS ts = (TS) P_data; 787cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 797cb94ee6SHong Zhang PC pc = cvode->pc; 807cb94ee6SHong Zhang Vec rr = cvode->w1,zz = cvode->w2; 817cb94ee6SHong Zhang PetscErrorCode ierr; 827cb94ee6SHong Zhang PetscScalar *r_data,*z_data; 837cb94ee6SHong Zhang 847cb94ee6SHong Zhang PetscFunctionBegin; 857cb94ee6SHong Zhang /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/ 867cb94ee6SHong Zhang r_data = (PetscScalar *) N_VGetArrayPointer(r); 877cb94ee6SHong Zhang z_data = (PetscScalar *) N_VGetArrayPointer(z); 887cb94ee6SHong Zhang ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr); 897cb94ee6SHong Zhang ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr); 904ac7836bSHong Zhang 917cb94ee6SHong Zhang /* Solve the Px=r and put the result in zz */ 927cb94ee6SHong Zhang ierr = PCApply(pc,rr,zz); CHKERRQ(ierr); 937cb94ee6SHong Zhang ierr = VecResetArray(rr); CHKERRQ(ierr); 947cb94ee6SHong Zhang ierr = VecResetArray(zz); CHKERRQ(ierr); 957cb94ee6SHong Zhang PetscFunctionReturn(0); 967cb94ee6SHong Zhang } 977cb94ee6SHong Zhang 987cb94ee6SHong Zhang /* 997cb94ee6SHong Zhang TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side. 1007cb94ee6SHong Zhang */ 1017cb94ee6SHong Zhang #undef __FUNCT__ 1027cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials" 1034ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx) 1047cb94ee6SHong Zhang { 1057cb94ee6SHong Zhang TS ts = (TS) ctx; 1067adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)ts)->comm; 1077cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1087cb94ee6SHong Zhang Vec yy = cvode->w1,yyd = cvode->w2; 1097cb94ee6SHong Zhang PetscScalar *y_data,*ydot_data; 1107cb94ee6SHong Zhang PetscErrorCode ierr; 1117cb94ee6SHong Zhang 1127cb94ee6SHong Zhang PetscFunctionBegin; 1135ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 1147cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(y); 1157cb94ee6SHong Zhang ydot_data = (PetscScalar *) N_VGetArrayPointer(ydot); 1164ac7836bSHong Zhang ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr); 1174ac7836bSHong Zhang ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr); 1184ac7836bSHong Zhang 1197cb94ee6SHong Zhang /* now compute the right hand side function */ 1207cb94ee6SHong Zhang ierr = TSComputeRHSFunction(ts,t,yy,yyd); CHKERRABORT(comm,ierr); 1217cb94ee6SHong Zhang ierr = VecResetArray(yy); CHKERRABORT(comm,ierr); 1227cb94ee6SHong Zhang ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr); 1234ac7836bSHong Zhang PetscFunctionReturn(0); 1247cb94ee6SHong Zhang } 1257cb94ee6SHong Zhang 1267cb94ee6SHong Zhang /* 1277cb94ee6SHong Zhang TSStep_Sundials_Nonlinear - Calls Sundials to integrate the ODE. 1287cb94ee6SHong Zhang */ 1297cb94ee6SHong Zhang #undef __FUNCT__ 1307cb94ee6SHong Zhang #define __FUNCT__ "TSStep_Sundials_Nonlinear" 1317cb94ee6SHong Zhang /* 1327cb94ee6SHong Zhang TSStep_Sundials_Nonlinear - 1337cb94ee6SHong Zhang 1347cb94ee6SHong Zhang steps - number of time steps 1357cb94ee6SHong Zhang time - time that integrater is terminated. 1367cb94ee6SHong Zhang */ 1377cb94ee6SHong Zhang PetscErrorCode TSStep_Sundials_Nonlinear(TS ts,int *steps,double *time) 1387cb94ee6SHong Zhang { 1397cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1407cb94ee6SHong Zhang Vec sol = ts->vec_sol; 1417cb94ee6SHong Zhang PetscErrorCode ierr; 1425d47aee6SHong Zhang PetscInt i,max_steps = ts->max_steps,flag; 1437cb94ee6SHong Zhang long int its; 1447cb94ee6SHong Zhang realtype t,tout; 1457cb94ee6SHong Zhang PetscScalar *y_data; 1467cb94ee6SHong Zhang void *mem; 1475d47aee6SHong Zhang const PCType pctype; 1485d47aee6SHong Zhang PetscTruth pcnone; 1497cb94ee6SHong Zhang 1507cb94ee6SHong Zhang PetscFunctionBegin; 15128adb3f7SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 1527cb94ee6SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 15328adb3f7SHong Zhang if (!mem) SETERRQ(PETSC_ERR_MEM,"CVodeCreate() fails"); 1547cb94ee6SHong Zhang 15528adb3f7SHong Zhang /* Set the pointer to user-defined data */ 15628adb3f7SHong Zhang flag = CVodeSetUserData(mem, ts); 15728adb3f7SHong Zhang if (flag) SETERRQ(PETSC_ERR_ARG_BADPTR,"CVodeSetUserData() fails"); 15828adb3f7SHong Zhang 15928adb3f7SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 16028adb3f7SHong Zhang * user's right hand side function in u'=f(t,u), the inital time T0, and 16128adb3f7SHong Zhang * the initial dependent variable vector cvode->y */ 16228adb3f7SHong Zhang flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 16328adb3f7SHong Zhang if (flag){ 16428adb3f7SHong Zhang SETERRQ1(1,"CVodeInit() fails, flag %d",flag); 16528adb3f7SHong Zhang } 16628adb3f7SHong Zhang 16728adb3f7SHong Zhang flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 16828adb3f7SHong Zhang if (flag){ 16928adb3f7SHong Zhang SETERRQ1(1,"CVodeSStolerances() fails, flag %d",flag); 17028adb3f7SHong Zhang } 17128adb3f7SHong Zhang /* add this call when there is a request 17228adb3f7SHong Zhang flag = CVodeSVtolerances(mem,cvode->reltol,cvode->abstol); 17328adb3f7SHong Zhang if (flag){ 17428adb3f7SHong Zhang SETERRQ1(1,"CVodeSVtolerances() fails, flag %d",flag); 17528adb3f7SHong Zhang } 1767cb94ee6SHong Zhang */ 1777cb94ee6SHong Zhang 1787cb94ee6SHong Zhang /* initialize the number of steps */ 1797cb94ee6SHong Zhang *steps = -ts->steps; 1807cb94ee6SHong Zhang ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1817cb94ee6SHong Zhang 1827cb94ee6SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 1837cb94ee6SHong Zhang /* setup the ode integrator with the given preconditioner */ 1845d47aee6SHong Zhang ierr = PCGetType(cvode->pc,&pctype);CHKERRQ(ierr); 1855d47aee6SHong Zhang ierr = PetscTypeCompare((PetscObject)cvode->pc,PCNONE,&pcnone);CHKERRQ(ierr); 1865d47aee6SHong Zhang if (pcnone){ 1875d47aee6SHong Zhang flag = CVSpgmr(mem,PREC_NONE,0); 1885d47aee6SHong Zhang } else { 1897cb94ee6SHong Zhang flag = CVSpgmr(mem,PREC_LEFT,0); 1905d47aee6SHong Zhang } 19128adb3f7SHong Zhang if (flag) SETERRQ1(1,"CVSpgmr() fails, flag %d",flag); 1924ac7836bSHong Zhang 1937cb94ee6SHong Zhang /* Set preconditioner setup and solve routines Precond and PSolve, 1947cb94ee6SHong Zhang and the pointer to the user-defined block data */ 19528adb3f7SHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 19628adb3f7SHong Zhang if (flag) SETERRQ1(1,"CVSpilsSetPreconditioner() fails, flag %d", flag); 1977cb94ee6SHong Zhang 1985d47aee6SHong Zhang flag = CVSpilsSetGSType(mem, MODIFIED_GS); 19928adb3f7SHong Zhang if (flag) SETERRQ1(1,"CVSpgmrSetGSType() fails, flag %d",flag); 2005d47aee6SHong Zhang 2017cb94ee6SHong Zhang tout = ts->max_time; 2027cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr); 2037cb94ee6SHong Zhang N_VSetArrayPointer((realtype *)y_data,cvode->y); 2047cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 2057cb94ee6SHong Zhang for (i = 0; i < max_steps; i++) { 2067cb94ee6SHong Zhang if (ts->ptime >= ts->max_time) break; 207*2bfc04deSHong Zhang if (cvode->monitorstep){ 20828adb3f7SHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP); 209*2bfc04deSHong Zhang } else { 210*2bfc04deSHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL); 211*2bfc04deSHong Zhang } 21228adb3f7SHong Zhang if (flag)SETERRQ1(1,"CVode() fails, flag %d",flag); 2137cb94ee6SHong Zhang if (t > ts->max_time && cvode->exact_final_time) { 2147cb94ee6SHong Zhang /* interpolate to final requested time */ 2157cb94ee6SHong Zhang ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr); 2167cb94ee6SHong Zhang t = tout; 2177cb94ee6SHong Zhang } 2187cb94ee6SHong Zhang ts->time_step = t - ts->ptime; 2197cb94ee6SHong Zhang ts->ptime = t; 2207cb94ee6SHong Zhang 2217cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 2227cb94ee6SHong Zhang ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr); 2237cb94ee6SHong Zhang ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr); 2247cb94ee6SHong Zhang ierr = VecResetArray(cvode->w1); CHKERRQ(ierr); 2257cb94ee6SHong Zhang ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr); 2267cb94ee6SHong Zhang ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr); 2277cb94ee6SHong Zhang ts->nonlinear_its = its; 2284ac7836bSHong Zhang ierr = CVSpilsGetNumLinIters(mem, &its); 2297cb94ee6SHong Zhang ts->linear_its = its; 2307cb94ee6SHong Zhang ts->steps++; 2317cb94ee6SHong Zhang ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr); 2327cb94ee6SHong Zhang } 2332c823083SHong Zhang cvode->mem = mem; 2347cb94ee6SHong Zhang *steps += ts->steps; 2357cb94ee6SHong Zhang *time = t; 2367cb94ee6SHong Zhang PetscFunctionReturn(0); 2377cb94ee6SHong Zhang } 2387cb94ee6SHong Zhang 2397cb94ee6SHong Zhang #undef __FUNCT__ 2407cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials" 2417cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts) 2427cb94ee6SHong Zhang { 2437cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2447cb94ee6SHong Zhang PetscErrorCode ierr; 2457cb94ee6SHong Zhang 2467cb94ee6SHong Zhang PetscFunctionBegin; 2477cb94ee6SHong Zhang if (cvode->pmat) {ierr = MatDestroy(cvode->pmat);CHKERRQ(ierr);} 2487cb94ee6SHong Zhang if (cvode->pc) {ierr = PCDestroy(cvode->pc);CHKERRQ(ierr);} 2497cb94ee6SHong Zhang if (cvode->update) {ierr = VecDestroy(cvode->update);CHKERRQ(ierr);} 2507cb94ee6SHong Zhang if (cvode->func) {ierr = VecDestroy(cvode->func);CHKERRQ(ierr);} 2517cb94ee6SHong Zhang if (cvode->rhs) {ierr = VecDestroy(cvode->rhs);CHKERRQ(ierr);} 2527cb94ee6SHong Zhang if (cvode->w1) {ierr = VecDestroy(cvode->w1);CHKERRQ(ierr);} 2537cb94ee6SHong Zhang if (cvode->w2) {ierr = VecDestroy(cvode->w2);CHKERRQ(ierr);} 254a07356b8SSatish Balay ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr); 2552c823083SHong Zhang if (cvode->mem) {CVodeFree(&cvode->mem);} 2567cb94ee6SHong Zhang ierr = PetscFree(cvode);CHKERRQ(ierr); 2577cb94ee6SHong Zhang PetscFunctionReturn(0); 2587cb94ee6SHong Zhang } 2597cb94ee6SHong Zhang 2607cb94ee6SHong Zhang #undef __FUNCT__ 2617cb94ee6SHong Zhang #define __FUNCT__ "TSSetUp_Sundials_Nonlinear" 2627cb94ee6SHong Zhang PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts) 2637cb94ee6SHong Zhang { 2647cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2657cb94ee6SHong Zhang PetscErrorCode ierr; 2667cb94ee6SHong Zhang int glosize,locsize,i; 2677cb94ee6SHong Zhang PetscScalar *y_data,*parray; 2687cb94ee6SHong Zhang 2697cb94ee6SHong Zhang PetscFunctionBegin; 2707cb94ee6SHong Zhang ierr = PCSetFromOptions(cvode->pc);CHKERRQ(ierr); 2717cb94ee6SHong Zhang /* get the vector size */ 2727cb94ee6SHong Zhang ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr); 2737cb94ee6SHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr); 2747cb94ee6SHong Zhang 2757cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 276a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 2777cb94ee6SHong Zhang if (!cvode->y) SETERRQ(1,"cvode->y is not allocated"); 2787cb94ee6SHong Zhang 27928adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 2807cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr); 2817cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y); 2827cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 2837cb94ee6SHong Zhang /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/ 2847cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 2857cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr); 2867cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->func);CHKERRQ(ierr); 2877cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr); 2887cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->func);CHKERRQ(ierr); 2897cb94ee6SHong Zhang 2907cb94ee6SHong Zhang /* 2917cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 2927cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 2937cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 2947cb94ee6SHong Zhang */ 2957adad957SLisandro Dalcin ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr); 2967adad957SLisandro Dalcin ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr); 2977cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr); 2987cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr); 2997cb94ee6SHong Zhang PetscFunctionReturn(0); 3007cb94ee6SHong Zhang } 3017cb94ee6SHong Zhang 3026fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 303dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0}; 3046fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 305dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0}; 306a04cf4d8SBarry Smith 3077cb94ee6SHong Zhang #undef __FUNCT__ 3087cb94ee6SHong Zhang #define __FUNCT__ "TSSetFromOptions_Sundials_Nonlinear" 3097cb94ee6SHong Zhang PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts) 3107cb94ee6SHong Zhang { 3117cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 3127cb94ee6SHong Zhang PetscErrorCode ierr; 3137cb94ee6SHong Zhang int indx; 3147cb94ee6SHong Zhang PetscTruth flag; 3157cb94ee6SHong Zhang 3167cb94ee6SHong Zhang PetscFunctionBegin; 3177cb94ee6SHong Zhang ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr); 3186fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr); 3197cb94ee6SHong Zhang if (flag) { 3206fadb2cdSHong Zhang ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr); 3217cb94ee6SHong Zhang } 3226fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr); 3237cb94ee6SHong Zhang if (flag) { 3247cb94ee6SHong Zhang ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 3257cb94ee6SHong Zhang } 3267cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr); 3277cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr); 3287cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr); 3297cb94ee6SHong Zhang ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr); 33090d69ab7SBarry Smith ierr = PetscOptionsTruth("-ts_sundials_exact_final_time","Allow SUNDIALS to stop near the final time, not exactly on it","TSSundialsSetExactFinalTime",cvode->exact_final_time,&cvode->exact_final_time,PETSC_NULL);CHKERRQ(ierr); 331*2bfc04deSHong Zhang ierr = PetscOptionsTruth("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr); 3327cb94ee6SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 3337cb94ee6SHong Zhang PetscFunctionReturn(0); 3347cb94ee6SHong Zhang } 3357cb94ee6SHong Zhang 3367cb94ee6SHong Zhang #undef __FUNCT__ 3377cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials" 3387cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 3397cb94ee6SHong Zhang { 3407cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 3417cb94ee6SHong Zhang PetscErrorCode ierr; 3427cb94ee6SHong Zhang char *type; 3437cb94ee6SHong Zhang char atype[] = "Adams"; 3447cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 3457cb94ee6SHong Zhang PetscTruth iascii,isstring; 3462c823083SHong Zhang long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 3472c823083SHong Zhang PetscInt qlast,qcur; 3485d47aee6SHong Zhang PetscReal hinused,hlast,hcur,tcur,tolsfac; 3497cb94ee6SHong Zhang 3507cb94ee6SHong Zhang PetscFunctionBegin; 3517cb94ee6SHong Zhang if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;} 3527cb94ee6SHong Zhang else {type = btype;} 3537cb94ee6SHong Zhang 3547cb94ee6SHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 3557cb94ee6SHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 3567cb94ee6SHong Zhang if (iascii) { 3577cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr); 3587cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr); 3597cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr); 3607cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr); 3617cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr); 3627cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 3637cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 3647cb94ee6SHong Zhang } else { 3657cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 3667cb94ee6SHong Zhang } 3672c823083SHong Zhang 3685d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 3695d47aee6SHong Zhang ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr); 3705d47aee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr); 3712c823083SHong Zhang ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals, 3722c823083SHong Zhang &nlinsetups,&nfails,&qlast,&qcur, 3732c823083SHong Zhang &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr); 3742c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr); 3752c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr); 3762c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr); 3772c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr); 3782c823083SHong Zhang 3792c823083SHong Zhang ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr); 3802c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr); 3812c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr); 3822c823083SHong Zhang 3832c823083SHong Zhang ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */ 3842c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr); 3852c823083SHong Zhang ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr); 3862c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr); 3872c823083SHong Zhang ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr); 3882c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr); 3892c823083SHong Zhang ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr); 3902c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr); 3912c823083SHong Zhang ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr); 3922c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr); 3932c823083SHong Zhang ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr); 3942c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr); 3957cb94ee6SHong Zhang } else if (isstring) { 3967cb94ee6SHong Zhang ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 3977cb94ee6SHong Zhang } else { 3987cb94ee6SHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name); 3997cb94ee6SHong Zhang } 4007cb94ee6SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 4017cb94ee6SHong Zhang ierr = PCView(cvode->pc,viewer);CHKERRQ(ierr); 4027cb94ee6SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 4037cb94ee6SHong Zhang PetscFunctionReturn(0); 4047cb94ee6SHong Zhang } 4057cb94ee6SHong Zhang 4067cb94ee6SHong Zhang 4077cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 4087cb94ee6SHong Zhang EXTERN_C_BEGIN 4097cb94ee6SHong Zhang #undef __FUNCT__ 4107cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials" 4116fadb2cdSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 4127cb94ee6SHong Zhang { 4137cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4147cb94ee6SHong Zhang 4157cb94ee6SHong Zhang PetscFunctionBegin; 4167cb94ee6SHong Zhang cvode->cvode_type = type; 4177cb94ee6SHong Zhang PetscFunctionReturn(0); 4187cb94ee6SHong Zhang } 4197cb94ee6SHong Zhang EXTERN_C_END 4207cb94ee6SHong Zhang 4217cb94ee6SHong Zhang EXTERN_C_BEGIN 4227cb94ee6SHong Zhang #undef __FUNCT__ 4237cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials" 4247cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart_Sundials(TS ts,int restart) 4257cb94ee6SHong Zhang { 4267cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4277cb94ee6SHong Zhang 4287cb94ee6SHong Zhang PetscFunctionBegin; 4297cb94ee6SHong Zhang cvode->restart = restart; 4307cb94ee6SHong Zhang PetscFunctionReturn(0); 4317cb94ee6SHong Zhang } 4327cb94ee6SHong Zhang EXTERN_C_END 4337cb94ee6SHong Zhang 4347cb94ee6SHong Zhang EXTERN_C_BEGIN 4357cb94ee6SHong Zhang #undef __FUNCT__ 4367cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials" 4377cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance_Sundials(TS ts,double tol) 4387cb94ee6SHong Zhang { 4397cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4407cb94ee6SHong Zhang 4417cb94ee6SHong Zhang PetscFunctionBegin; 4427cb94ee6SHong Zhang cvode->linear_tol = tol; 4437cb94ee6SHong Zhang PetscFunctionReturn(0); 4447cb94ee6SHong Zhang } 4457cb94ee6SHong Zhang EXTERN_C_END 4467cb94ee6SHong Zhang 4477cb94ee6SHong Zhang EXTERN_C_BEGIN 4487cb94ee6SHong Zhang #undef __FUNCT__ 4497cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials" 4507cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 4517cb94ee6SHong Zhang { 4527cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4537cb94ee6SHong Zhang 4547cb94ee6SHong Zhang PetscFunctionBegin; 4557cb94ee6SHong Zhang cvode->gtype = type; 4567cb94ee6SHong Zhang PetscFunctionReturn(0); 4577cb94ee6SHong Zhang } 4587cb94ee6SHong Zhang EXTERN_C_END 4597cb94ee6SHong Zhang 4607cb94ee6SHong Zhang EXTERN_C_BEGIN 4617cb94ee6SHong Zhang #undef __FUNCT__ 4627cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials" 4637cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel) 4647cb94ee6SHong Zhang { 4657cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4667cb94ee6SHong Zhang 4677cb94ee6SHong Zhang PetscFunctionBegin; 4687cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 4697cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 4707cb94ee6SHong Zhang PetscFunctionReturn(0); 4717cb94ee6SHong Zhang } 4727cb94ee6SHong Zhang EXTERN_C_END 4737cb94ee6SHong Zhang 4747cb94ee6SHong Zhang EXTERN_C_BEGIN 4757cb94ee6SHong Zhang #undef __FUNCT__ 4767cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials" 4777cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC_Sundials(TS ts,PC *pc) 4787cb94ee6SHong Zhang { 4797cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4807cb94ee6SHong Zhang 4817cb94ee6SHong Zhang PetscFunctionBegin; 4827cb94ee6SHong Zhang *pc = cvode->pc; 4837cb94ee6SHong Zhang PetscFunctionReturn(0); 4847cb94ee6SHong Zhang } 4857cb94ee6SHong Zhang EXTERN_C_END 4867cb94ee6SHong Zhang 4877cb94ee6SHong Zhang EXTERN_C_BEGIN 4887cb94ee6SHong Zhang #undef __FUNCT__ 4897cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials" 4907cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 4917cb94ee6SHong Zhang { 4927cb94ee6SHong Zhang PetscFunctionBegin; 4932c823083SHong Zhang if (nonlin) *nonlin = ts->nonlinear_its; 4942c823083SHong Zhang if (lin) *lin = ts->linear_its; 4957cb94ee6SHong Zhang PetscFunctionReturn(0); 4967cb94ee6SHong Zhang } 4977cb94ee6SHong Zhang EXTERN_C_END 4987cb94ee6SHong Zhang 4997cb94ee6SHong Zhang EXTERN_C_BEGIN 5007cb94ee6SHong Zhang #undef __FUNCT__ 5017cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials" 5027cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime_Sundials(TS ts,PetscTruth s) 5037cb94ee6SHong Zhang { 5047cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5057cb94ee6SHong Zhang 5067cb94ee6SHong Zhang PetscFunctionBegin; 5077cb94ee6SHong Zhang cvode->exact_final_time = s; 5087cb94ee6SHong Zhang PetscFunctionReturn(0); 5097cb94ee6SHong Zhang } 5107cb94ee6SHong Zhang EXTERN_C_END 511*2bfc04deSHong Zhang 512*2bfc04deSHong Zhang EXTERN_C_BEGIN 513*2bfc04deSHong Zhang #undef __FUNCT__ 514*2bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials" 515*2bfc04deSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscTruth s) 516*2bfc04deSHong Zhang { 517*2bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 518*2bfc04deSHong Zhang 519*2bfc04deSHong Zhang PetscFunctionBegin; 520*2bfc04deSHong Zhang cvode->monitorstep = s; 521*2bfc04deSHong Zhang PetscFunctionReturn(0); 522*2bfc04deSHong Zhang } 523*2bfc04deSHong Zhang EXTERN_C_END 5247cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 5257cb94ee6SHong Zhang 5267cb94ee6SHong Zhang #undef __FUNCT__ 5277cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations" 5287cb94ee6SHong Zhang /*@C 5297cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 5307cb94ee6SHong Zhang 5317cb94ee6SHong Zhang Not Collective 5327cb94ee6SHong Zhang 5337cb94ee6SHong Zhang Input parameters: 5347cb94ee6SHong Zhang . ts - the time-step context 5357cb94ee6SHong Zhang 5367cb94ee6SHong Zhang Output Parameters: 5377cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 5387cb94ee6SHong Zhang - lin - number of linear iterations 5397cb94ee6SHong Zhang 5407cb94ee6SHong Zhang Level: advanced 5417cb94ee6SHong Zhang 5427cb94ee6SHong Zhang Notes: 5437cb94ee6SHong Zhang These return the number since the creation of the TS object 5447cb94ee6SHong Zhang 5457cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations 5467cb94ee6SHong Zhang 5477cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(), 5487cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 5497cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 5507cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 5517cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 5527cb94ee6SHong Zhang 5537cb94ee6SHong Zhang @*/ 5547cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 5557cb94ee6SHong Zhang { 5567cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,int*,int*); 5577cb94ee6SHong Zhang 5587cb94ee6SHong Zhang PetscFunctionBegin; 5597cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetIterations_C",(void (**)(void))&f);CHKERRQ(ierr); 5607cb94ee6SHong Zhang if (f) { 5617cb94ee6SHong Zhang ierr = (*f)(ts,nonlin,lin);CHKERRQ(ierr); 5627cb94ee6SHong Zhang } 5637cb94ee6SHong Zhang PetscFunctionReturn(0); 5647cb94ee6SHong Zhang } 5657cb94ee6SHong Zhang 5667cb94ee6SHong Zhang #undef __FUNCT__ 5677cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType" 5687cb94ee6SHong Zhang /*@ 5697cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 5707cb94ee6SHong Zhang 5717cb94ee6SHong Zhang Collective on TS 5727cb94ee6SHong Zhang 5737cb94ee6SHong Zhang Input parameters: 5747cb94ee6SHong Zhang + ts - the time-step context 5757cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 5767cb94ee6SHong Zhang 5777cb94ee6SHong Zhang Level: intermediate 5787cb94ee6SHong Zhang 5797cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula 5807cb94ee6SHong Zhang 5817cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetGMRESRestart(), 5827cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 5837cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 5847cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 5857cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 5867cb94ee6SHong Zhang @*/ 5876fadb2cdSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType(TS ts,TSSundialsLmmType type) 5887cb94ee6SHong Zhang { 5896fadb2cdSHong Zhang PetscErrorCode ierr,(*f)(TS,TSSundialsLmmType); 5907cb94ee6SHong Zhang 5917cb94ee6SHong Zhang PetscFunctionBegin; 5927cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 5937cb94ee6SHong Zhang if (f) { 5947cb94ee6SHong Zhang ierr = (*f)(ts,type);CHKERRQ(ierr); 5957cb94ee6SHong Zhang } 5967cb94ee6SHong Zhang PetscFunctionReturn(0); 5977cb94ee6SHong Zhang } 5987cb94ee6SHong Zhang 5997cb94ee6SHong Zhang #undef __FUNCT__ 6007cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart" 6017cb94ee6SHong Zhang /*@ 6027cb94ee6SHong Zhang TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by 6037cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 6047cb94ee6SHong Zhang this is ALSO the maximum number of GMRES steps that will be used. 6057cb94ee6SHong Zhang 6067cb94ee6SHong Zhang Collective on TS 6077cb94ee6SHong Zhang 6087cb94ee6SHong Zhang Input parameters: 6097cb94ee6SHong Zhang + ts - the time-step context 6107cb94ee6SHong Zhang - restart - number of direction vectors (the restart size). 6117cb94ee6SHong Zhang 6127cb94ee6SHong Zhang Level: advanced 6137cb94ee6SHong Zhang 6147cb94ee6SHong Zhang .keywords: GMRES, restart 6157cb94ee6SHong Zhang 6167cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 6177cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 6187cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6197cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 6207cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 6217cb94ee6SHong Zhang 6227cb94ee6SHong Zhang @*/ 6237cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart(TS ts,int restart) 6247cb94ee6SHong Zhang { 6257cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,int); 6267cb94ee6SHong Zhang 6277cb94ee6SHong Zhang PetscFunctionBegin; 6287cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGMRESRestart_C",(void (**)(void))&f);CHKERRQ(ierr); 6297cb94ee6SHong Zhang if (f) { 6307cb94ee6SHong Zhang ierr = (*f)(ts,restart);CHKERRQ(ierr); 6317cb94ee6SHong Zhang } 6327cb94ee6SHong Zhang PetscFunctionReturn(0); 6337cb94ee6SHong Zhang } 6347cb94ee6SHong Zhang 6357cb94ee6SHong Zhang #undef __FUNCT__ 6367cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance" 6377cb94ee6SHong Zhang /*@ 6387cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 6397cb94ee6SHong Zhang system by SUNDIALS. 6407cb94ee6SHong Zhang 6417cb94ee6SHong Zhang Collective on TS 6427cb94ee6SHong Zhang 6437cb94ee6SHong Zhang Input parameters: 6447cb94ee6SHong Zhang + ts - the time-step context 6457cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 6467cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 6477cb94ee6SHong Zhang 6487cb94ee6SHong Zhang Level: advanced 6497cb94ee6SHong Zhang 6507cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS 6517cb94ee6SHong Zhang 6527cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6537cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 6547cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6557cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 6567cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 6577cb94ee6SHong Zhang 6587cb94ee6SHong Zhang @*/ 6597cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance(TS ts,double tol) 6607cb94ee6SHong Zhang { 6617cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,double); 6627cb94ee6SHong Zhang 6637cb94ee6SHong Zhang PetscFunctionBegin; 6647cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 6657cb94ee6SHong Zhang if (f) { 6667cb94ee6SHong Zhang ierr = (*f)(ts,tol);CHKERRQ(ierr); 6677cb94ee6SHong Zhang } 6687cb94ee6SHong Zhang PetscFunctionReturn(0); 6697cb94ee6SHong Zhang } 6707cb94ee6SHong Zhang 6717cb94ee6SHong Zhang #undef __FUNCT__ 6727cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType" 6737cb94ee6SHong Zhang /*@ 6747cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 6757cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 6767cb94ee6SHong Zhang 6777cb94ee6SHong Zhang Collective on TS 6787cb94ee6SHong Zhang 6797cb94ee6SHong Zhang Input parameters: 6807cb94ee6SHong Zhang + ts - the time-step context 6817cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 6827cb94ee6SHong Zhang 6837cb94ee6SHong Zhang Level: advanced 6847cb94ee6SHong Zhang 6857cb94ee6SHong Zhang .keywords: Sundials, orthogonalization 6867cb94ee6SHong Zhang 6877cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6887cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 6897cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6907cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 6917cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 6927cb94ee6SHong Zhang 6937cb94ee6SHong Zhang @*/ 6947cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 6957cb94ee6SHong Zhang { 6967cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,TSSundialsGramSchmidtType); 6977cb94ee6SHong Zhang 6987cb94ee6SHong Zhang PetscFunctionBegin; 6997cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",(void (**)(void))&f);CHKERRQ(ierr); 7007cb94ee6SHong Zhang if (f) { 7017cb94ee6SHong Zhang ierr = (*f)(ts,type);CHKERRQ(ierr); 7027cb94ee6SHong Zhang } 7037cb94ee6SHong Zhang PetscFunctionReturn(0); 7047cb94ee6SHong Zhang } 7057cb94ee6SHong Zhang 7067cb94ee6SHong Zhang #undef __FUNCT__ 7077cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance" 7087cb94ee6SHong Zhang /*@ 7097cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 7107cb94ee6SHong Zhang Sundials for error control. 7117cb94ee6SHong Zhang 7127cb94ee6SHong Zhang Collective on TS 7137cb94ee6SHong Zhang 7147cb94ee6SHong Zhang Input parameters: 7157cb94ee6SHong Zhang + ts - the time-step context 7167cb94ee6SHong Zhang . aabs - the absolute tolerance 7177cb94ee6SHong Zhang - rel - the relative tolerance 7187cb94ee6SHong Zhang 7197cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 7207cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 7217cb94ee6SHong Zhang 7227cb94ee6SHong Zhang Level: intermediate 7237cb94ee6SHong Zhang 7247cb94ee6SHong Zhang .keywords: Sundials, tolerance 7257cb94ee6SHong Zhang 7267cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7277cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 7287cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7297cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 7307cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 7317cb94ee6SHong Zhang 7327cb94ee6SHong Zhang @*/ 7337cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance(TS ts,double aabs,double rel) 7347cb94ee6SHong Zhang { 7357cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,double,double); 7367cb94ee6SHong Zhang 7377cb94ee6SHong Zhang PetscFunctionBegin; 7387cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 7397cb94ee6SHong Zhang if (f) { 7407cb94ee6SHong Zhang ierr = (*f)(ts,aabs,rel);CHKERRQ(ierr); 7417cb94ee6SHong Zhang } 7427cb94ee6SHong Zhang PetscFunctionReturn(0); 7437cb94ee6SHong Zhang } 7447cb94ee6SHong Zhang 7457cb94ee6SHong Zhang #undef __FUNCT__ 7467cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC" 7477cb94ee6SHong Zhang /*@ 7487cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 7497cb94ee6SHong Zhang 7507cb94ee6SHong Zhang Input Parameter: 7517cb94ee6SHong Zhang . ts - the time-step context 7527cb94ee6SHong Zhang 7537cb94ee6SHong Zhang Output Parameter: 7547cb94ee6SHong Zhang . pc - the preconditioner context 7557cb94ee6SHong Zhang 7567cb94ee6SHong Zhang Level: advanced 7577cb94ee6SHong Zhang 7587cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7597cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 7607cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7617cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 7627cb94ee6SHong Zhang @*/ 7637cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC(TS ts,PC *pc) 7647cb94ee6SHong Zhang { 7657cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,PC *); 7667cb94ee6SHong Zhang 7677cb94ee6SHong Zhang PetscFunctionBegin; 7687cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetPC_C",(void (**)(void))&f);CHKERRQ(ierr); 7697cb94ee6SHong Zhang if (f) { 7707cb94ee6SHong Zhang ierr = (*f)(ts,pc);CHKERRQ(ierr); 7717cb94ee6SHong Zhang } else { 7727cb94ee6SHong Zhang SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"TS must be of Sundials type to extract the PC"); 7737cb94ee6SHong Zhang } 7747cb94ee6SHong Zhang 7757cb94ee6SHong Zhang PetscFunctionReturn(0); 7767cb94ee6SHong Zhang } 7777cb94ee6SHong Zhang 7787cb94ee6SHong Zhang #undef __FUNCT__ 7797cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime" 7807cb94ee6SHong Zhang /*@ 7817cb94ee6SHong Zhang TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the 7827cb94ee6SHong Zhang exact final time requested by the user or just returns it at the final time 7837cb94ee6SHong Zhang it computed. (Defaults to true). 7847cb94ee6SHong Zhang 7857cb94ee6SHong Zhang Input Parameter: 7867cb94ee6SHong Zhang + ts - the time-step context 7877cb94ee6SHong Zhang - ft - PETSC_TRUE if interpolates, else PETSC_FALSE 7887cb94ee6SHong Zhang 7897cb94ee6SHong Zhang Level: beginner 7907cb94ee6SHong Zhang 7917cb94ee6SHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7927cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 7937cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7947cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 7957cb94ee6SHong Zhang @*/ 7967cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime(TS ts,PetscTruth ft) 7977cb94ee6SHong Zhang { 7987cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,PetscTruth); 7997cb94ee6SHong Zhang 8007cb94ee6SHong Zhang PetscFunctionBegin; 8017cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetExactFinalTime_C",(void (**)(void))&f);CHKERRQ(ierr); 8027cb94ee6SHong Zhang if (f) { 8037cb94ee6SHong Zhang ierr = (*f)(ts,ft);CHKERRQ(ierr); 8047cb94ee6SHong Zhang } 8057cb94ee6SHong Zhang PetscFunctionReturn(0); 8067cb94ee6SHong Zhang } 8077cb94ee6SHong Zhang 808*2bfc04deSHong Zhang #undef __FUNCT__ 809*2bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps" 810*2bfc04deSHong Zhang /*@ 811*2bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 812*2bfc04deSHong Zhang 813*2bfc04deSHong Zhang Input Parameter: 814*2bfc04deSHong Zhang + ts - the time-step context 815*2bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 816*2bfc04deSHong Zhang 817*2bfc04deSHong Zhang Level: beginner 818*2bfc04deSHong Zhang 819*2bfc04deSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 820*2bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 821*2bfc04deSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 822*2bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 823*2bfc04deSHong Zhang @*/ 824*2bfc04deSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsMonitorInternalSteps(TS ts,PetscTruth ft) 825*2bfc04deSHong Zhang { 826*2bfc04deSHong Zhang PetscErrorCode ierr,(*f)(TS,PetscTruth); 827*2bfc04deSHong Zhang 828*2bfc04deSHong Zhang PetscFunctionBegin; 829*2bfc04deSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",(void (**)(void))&f);CHKERRQ(ierr); 830*2bfc04deSHong Zhang if (f) { 831*2bfc04deSHong Zhang ierr = (*f)(ts,ft);CHKERRQ(ierr); 832*2bfc04deSHong Zhang } 833*2bfc04deSHong Zhang PetscFunctionReturn(0); 834*2bfc04deSHong Zhang } 8357cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 8367cb94ee6SHong Zhang /*MC 8377cb94ee6SHong Zhang TS_Sundials - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 8387cb94ee6SHong Zhang 8397cb94ee6SHong Zhang Options Database: 8407cb94ee6SHong Zhang + -ts_sundials_type <bdf,adams> 8417cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 8427cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 8437cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 8447cb94ee6SHong Zhang . -ts_sundials_linear_tolerance <tol> 8457cb94ee6SHong Zhang . -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions 8467cb94ee6SHong Zhang - -ts_sundials_not_exact_final_time -Allow SUNDIALS to stop near the final time, not exactly on it 8477cb94ee6SHong Zhang 8487cb94ee6SHong Zhang Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply 8497cb94ee6SHong Zhang only PETSc PC options 8507cb94ee6SHong Zhang 8517cb94ee6SHong Zhang Level: beginner 8527cb94ee6SHong Zhang 8537cb94ee6SHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(), 8547cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime() 8557cb94ee6SHong Zhang 8567cb94ee6SHong Zhang M*/ 8577cb94ee6SHong Zhang EXTERN_C_BEGIN 8587cb94ee6SHong Zhang #undef __FUNCT__ 8597cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials" 8607cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Sundials(TS ts) 8617cb94ee6SHong Zhang { 8627cb94ee6SHong Zhang TS_Sundials *cvode; 8637cb94ee6SHong Zhang PetscErrorCode ierr; 8647cb94ee6SHong Zhang 8657cb94ee6SHong Zhang PetscFunctionBegin; 8667cb94ee6SHong Zhang if (ts->problem_type != TS_NONLINEAR) { 8677cb94ee6SHong Zhang SETERRQ(PETSC_ERR_SUP,"Only support for nonlinear problems"); 8687cb94ee6SHong Zhang } 86928adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 87028adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 8717cb94ee6SHong Zhang ts->ops->setup = TSSetUp_Sundials_Nonlinear; 8727cb94ee6SHong Zhang ts->ops->step = TSStep_Sundials_Nonlinear; 8737cb94ee6SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_Sundials_Nonlinear; 8747cb94ee6SHong Zhang 87538f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr); 8767adad957SLisandro Dalcin ierr = PCCreate(((PetscObject)ts)->comm,&cvode->pc);CHKERRQ(ierr); 8777cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->pc);CHKERRQ(ierr); 8787cb94ee6SHong Zhang ts->data = (void*)cvode; 8796fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 8806fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 8817cb94ee6SHong Zhang cvode->restart = 5; 8827cb94ee6SHong Zhang cvode->linear_tol = .05; 8837cb94ee6SHong Zhang 884*2bfc04deSHong Zhang cvode->exact_final_time = PETSC_TRUE; 885*2bfc04deSHong Zhang cvode->monitorstep = PETSC_FALSE; 8867cb94ee6SHong Zhang 887a07356b8SSatish Balay ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr); 8887cb94ee6SHong Zhang /* set tolerance for Sundials */ 8897cb94ee6SHong Zhang cvode->reltol = 1e-6; 8902c823083SHong Zhang cvode->abstol = 1e-6; 8917cb94ee6SHong Zhang 8927cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials", 8937cb94ee6SHong Zhang TSSundialsSetType_Sundials);CHKERRQ(ierr); 8947cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C", 8957cb94ee6SHong Zhang "TSSundialsSetGMRESRestart_Sundials", 8967cb94ee6SHong Zhang TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr); 8977cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C", 8987cb94ee6SHong Zhang "TSSundialsSetLinearTolerance_Sundials", 8997cb94ee6SHong Zhang TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 9007cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C", 9017cb94ee6SHong Zhang "TSSundialsSetGramSchmidtType_Sundials", 9027cb94ee6SHong Zhang TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 9037cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C", 9047cb94ee6SHong Zhang "TSSundialsSetTolerance_Sundials", 9057cb94ee6SHong Zhang TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 9067cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C", 9077cb94ee6SHong Zhang "TSSundialsGetPC_Sundials", 9087cb94ee6SHong Zhang TSSundialsGetPC_Sundials);CHKERRQ(ierr); 9097cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C", 9107cb94ee6SHong Zhang "TSSundialsGetIterations_Sundials", 9117cb94ee6SHong Zhang TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 9127cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C", 9137cb94ee6SHong Zhang "TSSundialsSetExactFinalTime_Sundials", 9147cb94ee6SHong Zhang TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr); 915*2bfc04deSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C", 916*2bfc04deSHong Zhang "TSSundialsMonitorInternalSteps_Sundials", 917*2bfc04deSHong Zhang TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr); 918*2bfc04deSHong Zhang 9197cb94ee6SHong Zhang PetscFunctionReturn(0); 9207cb94ee6SHong Zhang } 9217cb94ee6SHong Zhang EXTERN_C_END 9227cb94ee6SHong Zhang 9237cb94ee6SHong Zhang 9247cb94ee6SHong Zhang 9257cb94ee6SHong Zhang 9267cb94ee6SHong Zhang 9277cb94ee6SHong Zhang 9287cb94ee6SHong Zhang 9297cb94ee6SHong Zhang 9307cb94ee6SHong Zhang 9317cb94ee6SHong Zhang 932