17cb94ee6SHong Zhang #define PETSCTS_DLL 27cb94ee6SHong Zhang 37cb94ee6SHong Zhang /* 4*db268bfaSHong 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. 7*db268bfaSHong Zhang Reference: sundials-2.3.0/examples/cvode/parallel/cvkryx_p.c 87cb94ee6SHong Zhang */ 956a740aaSHong Zhang #include "sundials.h" /*I "petscts.h" I*/ 107cb94ee6SHong Zhang 117cb94ee6SHong Zhang /* 127cb94ee6SHong Zhang TSPrecond_Sundials - function that we provide to SUNDIALS to 137cb94ee6SHong Zhang evaluate the preconditioner. 147cb94ee6SHong Zhang */ 157cb94ee6SHong Zhang #undef __FUNCT__ 167cb94ee6SHong Zhang #define __FUNCT__ "TSPrecond_Sundials" 177cb94ee6SHong Zhang PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy, 187cb94ee6SHong Zhang booleantype jok,booleantype *jcurPtr, 197cb94ee6SHong Zhang realtype _gamma,void *P_data, 207cb94ee6SHong Zhang N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3) 217cb94ee6SHong Zhang { 227cb94ee6SHong Zhang TS ts = (TS) P_data; 237cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 247cb94ee6SHong Zhang PC pc = cvode->pc; 257cb94ee6SHong Zhang PetscErrorCode ierr; 267cb94ee6SHong Zhang Mat Jac = ts->B; 277cb94ee6SHong Zhang Vec yy = cvode->w1; 287cb94ee6SHong Zhang PetscScalar one = 1.0,gm; 297cb94ee6SHong Zhang MatStructure str = DIFFERENT_NONZERO_PATTERN; 307cb94ee6SHong Zhang PetscScalar *y_data; 317cb94ee6SHong Zhang 327cb94ee6SHong Zhang PetscFunctionBegin; 337cb94ee6SHong Zhang /* This allows us to construct preconditioners in-place if we like */ 347cb94ee6SHong Zhang ierr = MatSetUnfactored(Jac);CHKERRQ(ierr); 357cb94ee6SHong Zhang 367cb94ee6SHong Zhang /* jok - TRUE means reuse current Jacobian else recompute Jacobian */ 377cb94ee6SHong Zhang if (jok) { 387cb94ee6SHong Zhang ierr = MatCopy(cvode->pmat,Jac,str);CHKERRQ(ierr); 397cb94ee6SHong Zhang *jcurPtr = FALSE; 407cb94ee6SHong Zhang } else { 417cb94ee6SHong Zhang /* make PETSc vector yy point to SUNDIALS vector y */ 427cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(y); 437cb94ee6SHong Zhang ierr = VecPlaceArray(yy,y_data); CHKERRQ(ierr); 444ac7836bSHong Zhang 457cb94ee6SHong Zhang /* compute the Jacobian */ 467cb94ee6SHong Zhang ierr = TSComputeRHSJacobian(ts,ts->ptime,yy,&Jac,&Jac,&str);CHKERRQ(ierr); 477cb94ee6SHong Zhang ierr = VecResetArray(yy); CHKERRQ(ierr); 484ac7836bSHong Zhang 497cb94ee6SHong Zhang /* copy the Jacobian matrix */ 507cb94ee6SHong Zhang if (!cvode->pmat) { 517cb94ee6SHong Zhang ierr = MatDuplicate(Jac,MAT_COPY_VALUES,&cvode->pmat);CHKERRQ(ierr); 527cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->pmat);CHKERRQ(ierr); 537cb94ee6SHong Zhang } 547cb94ee6SHong 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 cvode->linear_solves++; 967cb94ee6SHong Zhang PetscFunctionReturn(0); 977cb94ee6SHong Zhang } 987cb94ee6SHong Zhang 997cb94ee6SHong Zhang /* 1007cb94ee6SHong Zhang TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side. 1017cb94ee6SHong Zhang */ 1027cb94ee6SHong Zhang #undef __FUNCT__ 1037cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials" 1044ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx) 1057cb94ee6SHong Zhang { 1067cb94ee6SHong Zhang TS ts = (TS) ctx; 1077cb94ee6SHong Zhang MPI_Comm comm = ts->comm; 1087cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1097cb94ee6SHong Zhang Vec yy = cvode->w1,yyd = cvode->w2; 1107cb94ee6SHong Zhang PetscScalar *y_data,*ydot_data; 1117cb94ee6SHong Zhang PetscErrorCode ierr; 1127cb94ee6SHong Zhang 1137cb94ee6SHong Zhang PetscFunctionBegin; 1145ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 1157cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(y); 1167cb94ee6SHong Zhang ydot_data = (PetscScalar *) N_VGetArrayPointer(ydot); 1174ac7836bSHong Zhang ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr); 1184ac7836bSHong Zhang ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr); 1194ac7836bSHong Zhang 1207cb94ee6SHong Zhang /* now compute the right hand side function */ 1217cb94ee6SHong Zhang ierr = TSComputeRHSFunction(ts,t,yy,yyd); CHKERRABORT(comm,ierr); 1227cb94ee6SHong Zhang ierr = VecResetArray(yy); CHKERRABORT(comm,ierr); 1237cb94ee6SHong Zhang ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr); 1244ac7836bSHong Zhang PetscFunctionReturn(0); 1257cb94ee6SHong Zhang } 1267cb94ee6SHong Zhang 1277cb94ee6SHong Zhang /* 1287cb94ee6SHong Zhang TSStep_Sundials_Nonlinear - Calls Sundials to integrate the ODE. 1297cb94ee6SHong Zhang */ 1307cb94ee6SHong Zhang #undef __FUNCT__ 1317cb94ee6SHong Zhang #define __FUNCT__ "TSStep_Sundials_Nonlinear" 1327cb94ee6SHong Zhang /* 1337cb94ee6SHong Zhang TSStep_Sundials_Nonlinear - 1347cb94ee6SHong Zhang 1357cb94ee6SHong Zhang steps - number of time steps 1367cb94ee6SHong Zhang time - time that integrater is terminated. 1377cb94ee6SHong Zhang */ 1387cb94ee6SHong Zhang PetscErrorCode TSStep_Sundials_Nonlinear(TS ts,int *steps,double *time) 1397cb94ee6SHong Zhang { 1407cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1417cb94ee6SHong Zhang Vec sol = ts->vec_sol; 1427cb94ee6SHong Zhang PetscErrorCode ierr; 1437cb94ee6SHong Zhang int i,max_steps = ts->max_steps,flag; 1447cb94ee6SHong Zhang long int its; 1457cb94ee6SHong Zhang realtype t,tout; 1467cb94ee6SHong Zhang PetscScalar *y_data; 1477cb94ee6SHong Zhang void *mem; 1487cb94ee6SHong Zhang 1497cb94ee6SHong Zhang PetscFunctionBegin; 150*db268bfaSHong Zhang /* Call CVodeCreate to create the solver memory */ 1517cb94ee6SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 1527cb94ee6SHong Zhang if (!mem) SETERRQ(1,"CVodeCreate() fails"); 1537cb94ee6SHong Zhang flag = CVodeSetFdata(mem,ts); 1547cb94ee6SHong Zhang if (flag) SETERRQ(1,"CVodeSetFdata() fails"); 1557cb94ee6SHong Zhang 1567cb94ee6SHong Zhang /* 1577cb94ee6SHong Zhang Call CVodeMalloc to initialize the integrator memory: 1587cb94ee6SHong Zhang mem is the pointer to the integrator memory returned by CVodeCreate 1597cb94ee6SHong Zhang f is the user's right hand side function in y'=f(t,y) 1607cb94ee6SHong Zhang T0 is the initial time 1617cb94ee6SHong Zhang u is the initial dependent variable vector 1627cb94ee6SHong Zhang CV_SS specifies scalar relative and absolute tolerances 1637cb94ee6SHong Zhang reltol is the relative tolerance 1647cb94ee6SHong Zhang &abstol is a pointer to the scalar absolute tolerance 1657cb94ee6SHong Zhang */ 1667cb94ee6SHong Zhang flag = CVodeMalloc(mem,TSFunction_Sundials,ts->ptime,cvode->y,CV_SS,cvode->reltol,&cvode->abstol); 1677cb94ee6SHong Zhang if (flag) SETERRQ(1,"CVodeMalloc() fails"); 1687cb94ee6SHong Zhang 1697cb94ee6SHong Zhang /* initialize the number of steps */ 1707cb94ee6SHong Zhang *steps = -ts->steps; 1717cb94ee6SHong Zhang ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1727cb94ee6SHong Zhang 1737cb94ee6SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 1747cb94ee6SHong Zhang /* setup the ode integrator with the given preconditioner */ 1757cb94ee6SHong Zhang flag = CVSpgmr(mem,PREC_LEFT,0); 1767cb94ee6SHong Zhang if (flag) SETERRQ(1,"CVSpgmr() fails"); 1774ac7836bSHong Zhang 1784ac7836bSHong Zhang flag = CVSpilsSetGSType(mem, MODIFIED_GS); 1797cb94ee6SHong Zhang if (flag) SETERRQ(1,"CVSpgmrSetGSType() fails"); 1807cb94ee6SHong Zhang 1817cb94ee6SHong Zhang /* Set preconditioner setup and solve routines Precond and PSolve, 1827cb94ee6SHong Zhang and the pointer to the user-defined block data */ 1834ac7836bSHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials,ts); 1847cb94ee6SHong Zhang if (flag) SETERRQ(1,"CVSpgmrSetPreconditioner() fails"); 1857cb94ee6SHong Zhang 1867cb94ee6SHong Zhang tout = ts->max_time; 1877cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr); 1887cb94ee6SHong Zhang N_VSetArrayPointer((realtype *)y_data,cvode->y); 1897cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 1907cb94ee6SHong Zhang for (i = 0; i < max_steps; i++) { 1917cb94ee6SHong Zhang if (ts->ptime >= ts->max_time) break; 1927cb94ee6SHong Zhang ierr = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);CHKERRQ(ierr); 1937cb94ee6SHong Zhang ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr); 1947cb94ee6SHong Zhang cvode->nonlinear_solves += its; 1957cb94ee6SHong Zhang 1967cb94ee6SHong Zhang if (t > ts->max_time && cvode->exact_final_time) { 1977cb94ee6SHong Zhang /* interpolate to final requested time */ 1987cb94ee6SHong Zhang ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr); 1997cb94ee6SHong Zhang t = tout; 2007cb94ee6SHong Zhang } 2017cb94ee6SHong Zhang ts->time_step = t - ts->ptime; 2027cb94ee6SHong Zhang ts->ptime = t; 2037cb94ee6SHong Zhang 2047cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 2057cb94ee6SHong Zhang ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr); 2067cb94ee6SHong Zhang ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr); 2077cb94ee6SHong Zhang ierr = VecResetArray(cvode->w1); CHKERRQ(ierr); 2087cb94ee6SHong Zhang ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr); 2097cb94ee6SHong Zhang ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr); 2107cb94ee6SHong Zhang ts->nonlinear_its = its; 2114ac7836bSHong Zhang ierr = CVSpilsGetNumLinIters(mem, &its); 2127cb94ee6SHong Zhang ts->linear_its = its; 2137cb94ee6SHong Zhang ts->steps++; 2147cb94ee6SHong Zhang ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr); 2157cb94ee6SHong Zhang } 216f1945f72SSatish Balay CVodeFree(&mem); 2177cb94ee6SHong Zhang *steps += ts->steps; 2187cb94ee6SHong Zhang *time = t; 2197cb94ee6SHong Zhang PetscFunctionReturn(0); 2207cb94ee6SHong Zhang } 2217cb94ee6SHong Zhang 2227cb94ee6SHong Zhang #undef __FUNCT__ 2237cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials" 2247cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts) 2257cb94ee6SHong Zhang { 2267cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2277cb94ee6SHong Zhang PetscErrorCode ierr; 2287cb94ee6SHong Zhang 2297cb94ee6SHong Zhang PetscFunctionBegin; 2307cb94ee6SHong Zhang if (cvode->pmat) {ierr = MatDestroy(cvode->pmat);CHKERRQ(ierr);} 2317cb94ee6SHong Zhang if (cvode->pc) {ierr = PCDestroy(cvode->pc);CHKERRQ(ierr);} 2327cb94ee6SHong Zhang if (cvode->update) {ierr = VecDestroy(cvode->update);CHKERRQ(ierr);} 2337cb94ee6SHong Zhang if (cvode->func) {ierr = VecDestroy(cvode->func);CHKERRQ(ierr);} 2347cb94ee6SHong Zhang if (cvode->rhs) {ierr = VecDestroy(cvode->rhs);CHKERRQ(ierr);} 2357cb94ee6SHong Zhang if (cvode->w1) {ierr = VecDestroy(cvode->w1);CHKERRQ(ierr);} 2367cb94ee6SHong Zhang if (cvode->w2) {ierr = VecDestroy(cvode->w2);CHKERRQ(ierr);} 2377cb94ee6SHong Zhang ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr); 2387cb94ee6SHong Zhang ierr = PetscFree(cvode);CHKERRQ(ierr); 2397cb94ee6SHong Zhang PetscFunctionReturn(0); 2407cb94ee6SHong Zhang } 2417cb94ee6SHong Zhang 2427cb94ee6SHong Zhang #undef __FUNCT__ 2437cb94ee6SHong Zhang #define __FUNCT__ "TSSetUp_Sundials_Nonlinear" 2447cb94ee6SHong Zhang PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts) 2457cb94ee6SHong Zhang { 2467cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2477cb94ee6SHong Zhang PetscErrorCode ierr; 2487cb94ee6SHong Zhang int glosize,locsize,i; 2497cb94ee6SHong Zhang PetscScalar *y_data,*parray; 2507cb94ee6SHong Zhang 2517cb94ee6SHong Zhang PetscFunctionBegin; 2527cb94ee6SHong Zhang ierr = PCSetFromOptions(cvode->pc);CHKERRQ(ierr); 2537cb94ee6SHong Zhang /* get the vector size */ 2547cb94ee6SHong Zhang ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr); 2557cb94ee6SHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr); 2567cb94ee6SHong Zhang 2577cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 2587cb94ee6SHong Zhang cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 2597cb94ee6SHong Zhang if (!cvode->y) SETERRQ(1,"cvode->y is not allocated"); 2607cb94ee6SHong Zhang 2617cb94ee6SHong Zhang /* initialize N_Vec y */ 2627cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr); 2637cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y); 2647cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 2657cb94ee6SHong Zhang /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/ 2667cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 2677cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr); 2687cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->func);CHKERRQ(ierr); 2697cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr); 2707cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->func);CHKERRQ(ierr); 2717cb94ee6SHong Zhang 2727cb94ee6SHong Zhang /* 2737cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 2747cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 2757cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 2767cb94ee6SHong Zhang */ 2777cb94ee6SHong Zhang ierr = VecCreateMPIWithArray(ts->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr); 2787cb94ee6SHong Zhang ierr = VecCreateMPIWithArray(ts->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr); 2797cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr); 2807cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr); 2817cb94ee6SHong Zhang PetscFunctionReturn(0); 2827cb94ee6SHong Zhang } 2837cb94ee6SHong Zhang 2847cb94ee6SHong Zhang #undef __FUNCT__ 2857cb94ee6SHong Zhang #define __FUNCT__ "TSSetFromOptions_Sundials_Nonlinear" 2867cb94ee6SHong Zhang PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts) 2877cb94ee6SHong Zhang { 2887cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2897cb94ee6SHong Zhang PetscErrorCode ierr; 2907cb94ee6SHong Zhang int indx; 291*db268bfaSHong Zhang const char *btype[] = {" ","adams","bdf"},*otype[] = {"modified","unmodified"}; 2927cb94ee6SHong Zhang PetscTruth flag; 2937cb94ee6SHong Zhang 2947cb94ee6SHong Zhang PetscFunctionBegin; 2957cb94ee6SHong Zhang ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr); 296*db268bfaSHong Zhang ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",btype,3,"bdf",&indx,&flag);CHKERRQ(ierr); 2977cb94ee6SHong Zhang if (flag) { 2987cb94ee6SHong Zhang ierr = TSSundialsSetType(ts,(TSSundialsType)indx);CHKERRQ(ierr); 2997cb94ee6SHong Zhang } 3007cb94ee6SHong Zhang ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",otype,2,"unmodified",&indx,&flag);CHKERRQ(ierr); 3017cb94ee6SHong Zhang if (flag) { 3027cb94ee6SHong Zhang ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 3037cb94ee6SHong Zhang } 3047cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr); 3057cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr); 3067cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr); 3077cb94ee6SHong Zhang ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr); 308d7fc3247SBarry Smith ierr = PetscOptionsName("-ts_sundials_exact_final_time","Allow SUNDIALS to stop near the final time, not exactly on it","TSSundialsSetExactFinalTime",&flag);CHKERRQ(ierr); 309d7fc3247SBarry Smith if (flag) cvode->exact_final_time = PETSC_TRUE; 3107cb94ee6SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 3117cb94ee6SHong Zhang PetscFunctionReturn(0); 3127cb94ee6SHong Zhang } 3137cb94ee6SHong Zhang 3147cb94ee6SHong Zhang #undef __FUNCT__ 3157cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials" 3167cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 3177cb94ee6SHong Zhang { 3187cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 3197cb94ee6SHong Zhang PetscErrorCode ierr; 3207cb94ee6SHong Zhang char *type; 3217cb94ee6SHong Zhang char atype[] = "Adams"; 3227cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 3237cb94ee6SHong Zhang PetscTruth iascii,isstring; 3247cb94ee6SHong Zhang 3257cb94ee6SHong Zhang PetscFunctionBegin; 3267cb94ee6SHong Zhang if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;} 3277cb94ee6SHong Zhang else {type = btype;} 3287cb94ee6SHong Zhang 3297cb94ee6SHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 3307cb94ee6SHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 3317cb94ee6SHong Zhang if (iascii) { 3327cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr); 3337cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr); 3347cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr); 3357cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr); 3367cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr); 3377cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 3387cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 3397cb94ee6SHong Zhang } else { 3407cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 3417cb94ee6SHong Zhang } 3427cb94ee6SHong Zhang } else if (isstring) { 3437cb94ee6SHong Zhang ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 3447cb94ee6SHong Zhang } else { 3457cb94ee6SHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name); 3467cb94ee6SHong Zhang } 3477cb94ee6SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 3487cb94ee6SHong Zhang ierr = PCView(cvode->pc,viewer);CHKERRQ(ierr); 3497cb94ee6SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 3507cb94ee6SHong Zhang PetscFunctionReturn(0); 3517cb94ee6SHong Zhang } 3527cb94ee6SHong Zhang 3537cb94ee6SHong Zhang 3547cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 3557cb94ee6SHong Zhang EXTERN_C_BEGIN 3567cb94ee6SHong Zhang #undef __FUNCT__ 3577cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials" 3587cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType_Sundials(TS ts,TSSundialsType type) 3597cb94ee6SHong Zhang { 3607cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 3617cb94ee6SHong Zhang 3627cb94ee6SHong Zhang PetscFunctionBegin; 3637cb94ee6SHong Zhang cvode->cvode_type = type; 3647cb94ee6SHong Zhang PetscFunctionReturn(0); 3657cb94ee6SHong Zhang } 3667cb94ee6SHong Zhang EXTERN_C_END 3677cb94ee6SHong Zhang 3687cb94ee6SHong Zhang EXTERN_C_BEGIN 3697cb94ee6SHong Zhang #undef __FUNCT__ 3707cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials" 3717cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart_Sundials(TS ts,int restart) 3727cb94ee6SHong Zhang { 3737cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 3747cb94ee6SHong Zhang 3757cb94ee6SHong Zhang PetscFunctionBegin; 3767cb94ee6SHong Zhang cvode->restart = restart; 3777cb94ee6SHong Zhang PetscFunctionReturn(0); 3787cb94ee6SHong Zhang } 3797cb94ee6SHong Zhang EXTERN_C_END 3807cb94ee6SHong Zhang 3817cb94ee6SHong Zhang EXTERN_C_BEGIN 3827cb94ee6SHong Zhang #undef __FUNCT__ 3837cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials" 3847cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance_Sundials(TS ts,double tol) 3857cb94ee6SHong Zhang { 3867cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 3877cb94ee6SHong Zhang 3887cb94ee6SHong Zhang PetscFunctionBegin; 3897cb94ee6SHong Zhang cvode->linear_tol = tol; 3907cb94ee6SHong Zhang PetscFunctionReturn(0); 3917cb94ee6SHong Zhang } 3927cb94ee6SHong Zhang EXTERN_C_END 3937cb94ee6SHong Zhang 3947cb94ee6SHong Zhang EXTERN_C_BEGIN 3957cb94ee6SHong Zhang #undef __FUNCT__ 3967cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials" 3977cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 3987cb94ee6SHong Zhang { 3997cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4007cb94ee6SHong Zhang 4017cb94ee6SHong Zhang PetscFunctionBegin; 4027cb94ee6SHong Zhang cvode->gtype = type; 4037cb94ee6SHong Zhang PetscFunctionReturn(0); 4047cb94ee6SHong Zhang } 4057cb94ee6SHong Zhang EXTERN_C_END 4067cb94ee6SHong Zhang 4077cb94ee6SHong Zhang EXTERN_C_BEGIN 4087cb94ee6SHong Zhang #undef __FUNCT__ 4097cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials" 4107cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel) 4117cb94ee6SHong Zhang { 4127cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4137cb94ee6SHong Zhang 4147cb94ee6SHong Zhang PetscFunctionBegin; 4157cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 4167cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 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__ "TSSundialsGetPC_Sundials" 4247cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC_Sundials(TS ts,PC *pc) 4257cb94ee6SHong Zhang { 4267cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4277cb94ee6SHong Zhang 4287cb94ee6SHong Zhang PetscFunctionBegin; 4297cb94ee6SHong Zhang *pc = cvode->pc; 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__ "TSSundialsGetIterations_Sundials" 4377cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 4387cb94ee6SHong Zhang { 4397cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4407cb94ee6SHong Zhang 4417cb94ee6SHong Zhang PetscFunctionBegin; 4427cb94ee6SHong Zhang if (nonlin) *nonlin = cvode->nonlinear_solves; 4437cb94ee6SHong Zhang if (lin) *lin = cvode->linear_solves; 4447cb94ee6SHong Zhang PetscFunctionReturn(0); 4457cb94ee6SHong Zhang } 4467cb94ee6SHong Zhang EXTERN_C_END 4477cb94ee6SHong Zhang 4487cb94ee6SHong Zhang EXTERN_C_BEGIN 4497cb94ee6SHong Zhang #undef __FUNCT__ 4507cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials" 4517cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime_Sundials(TS ts,PetscTruth s) 4527cb94ee6SHong Zhang { 4537cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4547cb94ee6SHong Zhang 4557cb94ee6SHong Zhang PetscFunctionBegin; 4567cb94ee6SHong Zhang cvode->exact_final_time = s; 4577cb94ee6SHong Zhang PetscFunctionReturn(0); 4587cb94ee6SHong Zhang } 4597cb94ee6SHong Zhang EXTERN_C_END 4607cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 4617cb94ee6SHong Zhang 4627cb94ee6SHong Zhang #undef __FUNCT__ 4637cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations" 4647cb94ee6SHong Zhang /*@C 4657cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 4667cb94ee6SHong Zhang 4677cb94ee6SHong Zhang Not Collective 4687cb94ee6SHong Zhang 4697cb94ee6SHong Zhang Input parameters: 4707cb94ee6SHong Zhang . ts - the time-step context 4717cb94ee6SHong Zhang 4727cb94ee6SHong Zhang Output Parameters: 4737cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 4747cb94ee6SHong Zhang - lin - number of linear iterations 4757cb94ee6SHong Zhang 4767cb94ee6SHong Zhang Level: advanced 4777cb94ee6SHong Zhang 4787cb94ee6SHong Zhang Notes: 4797cb94ee6SHong Zhang These return the number since the creation of the TS object 4807cb94ee6SHong Zhang 4817cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations 4827cb94ee6SHong Zhang 4837cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(), 4847cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 4857cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 4867cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 4877cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 4887cb94ee6SHong Zhang 4897cb94ee6SHong Zhang @*/ 4907cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 4917cb94ee6SHong Zhang { 4927cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,int*,int*); 4937cb94ee6SHong Zhang 4947cb94ee6SHong Zhang PetscFunctionBegin; 4957cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetIterations_C",(void (**)(void))&f);CHKERRQ(ierr); 4967cb94ee6SHong Zhang if (f) { 4977cb94ee6SHong Zhang ierr = (*f)(ts,nonlin,lin);CHKERRQ(ierr); 4987cb94ee6SHong Zhang } 4997cb94ee6SHong Zhang PetscFunctionReturn(0); 5007cb94ee6SHong Zhang } 5017cb94ee6SHong Zhang 5027cb94ee6SHong Zhang #undef __FUNCT__ 5037cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType" 5047cb94ee6SHong Zhang /*@ 5057cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 5067cb94ee6SHong Zhang 5077cb94ee6SHong Zhang Collective on TS 5087cb94ee6SHong Zhang 5097cb94ee6SHong Zhang Input parameters: 5107cb94ee6SHong Zhang + ts - the time-step context 5117cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 5127cb94ee6SHong Zhang 5137cb94ee6SHong Zhang Level: intermediate 5147cb94ee6SHong Zhang 5157cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula 5167cb94ee6SHong Zhang 5177cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetGMRESRestart(), 5187cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 5197cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 5207cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 5217cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 5227cb94ee6SHong Zhang @*/ 5237cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType(TS ts,TSSundialsType type) 5247cb94ee6SHong Zhang { 5257cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,TSSundialsType); 5267cb94ee6SHong Zhang 5277cb94ee6SHong Zhang PetscFunctionBegin; 5287cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 5297cb94ee6SHong Zhang if (f) { 5307cb94ee6SHong Zhang ierr = (*f)(ts,type);CHKERRQ(ierr); 5317cb94ee6SHong Zhang } 5327cb94ee6SHong Zhang PetscFunctionReturn(0); 5337cb94ee6SHong Zhang } 5347cb94ee6SHong Zhang 5357cb94ee6SHong Zhang #undef __FUNCT__ 5367cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart" 5377cb94ee6SHong Zhang /*@ 5387cb94ee6SHong Zhang TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by 5397cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 5407cb94ee6SHong Zhang this is ALSO the maximum number of GMRES steps that will be used. 5417cb94ee6SHong Zhang 5427cb94ee6SHong Zhang Collective on TS 5437cb94ee6SHong Zhang 5447cb94ee6SHong Zhang Input parameters: 5457cb94ee6SHong Zhang + ts - the time-step context 5467cb94ee6SHong Zhang - restart - number of direction vectors (the restart size). 5477cb94ee6SHong Zhang 5487cb94ee6SHong Zhang Level: advanced 5497cb94ee6SHong Zhang 5507cb94ee6SHong Zhang .keywords: GMRES, restart 5517cb94ee6SHong Zhang 5527cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 5537cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 5547cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 5557cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 5567cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 5577cb94ee6SHong Zhang 5587cb94ee6SHong Zhang @*/ 5597cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart(TS ts,int restart) 5607cb94ee6SHong Zhang { 5617cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,int); 5627cb94ee6SHong Zhang 5637cb94ee6SHong Zhang PetscFunctionBegin; 5647cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGMRESRestart_C",(void (**)(void))&f);CHKERRQ(ierr); 5657cb94ee6SHong Zhang if (f) { 5667cb94ee6SHong Zhang ierr = (*f)(ts,restart);CHKERRQ(ierr); 5677cb94ee6SHong Zhang } 5687cb94ee6SHong Zhang PetscFunctionReturn(0); 5697cb94ee6SHong Zhang } 5707cb94ee6SHong Zhang 5717cb94ee6SHong Zhang #undef __FUNCT__ 5727cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance" 5737cb94ee6SHong Zhang /*@ 5747cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 5757cb94ee6SHong Zhang system by SUNDIALS. 5767cb94ee6SHong Zhang 5777cb94ee6SHong Zhang Collective on TS 5787cb94ee6SHong Zhang 5797cb94ee6SHong Zhang Input parameters: 5807cb94ee6SHong Zhang + ts - the time-step context 5817cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 5827cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 5837cb94ee6SHong Zhang 5847cb94ee6SHong Zhang Level: advanced 5857cb94ee6SHong Zhang 5867cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS 5877cb94ee6SHong Zhang 5887cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 5897cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 5907cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 5917cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 5927cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 5937cb94ee6SHong Zhang 5947cb94ee6SHong Zhang @*/ 5957cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance(TS ts,double tol) 5967cb94ee6SHong Zhang { 5977cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,double); 5987cb94ee6SHong Zhang 5997cb94ee6SHong Zhang PetscFunctionBegin; 6007cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 6017cb94ee6SHong Zhang if (f) { 6027cb94ee6SHong Zhang ierr = (*f)(ts,tol);CHKERRQ(ierr); 6037cb94ee6SHong Zhang } 6047cb94ee6SHong Zhang PetscFunctionReturn(0); 6057cb94ee6SHong Zhang } 6067cb94ee6SHong Zhang 6077cb94ee6SHong Zhang #undef __FUNCT__ 6087cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType" 6097cb94ee6SHong Zhang /*@ 6107cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 6117cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 6127cb94ee6SHong Zhang 6137cb94ee6SHong Zhang Collective on TS 6147cb94ee6SHong Zhang 6157cb94ee6SHong Zhang Input parameters: 6167cb94ee6SHong Zhang + ts - the time-step context 6177cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 6187cb94ee6SHong Zhang 6197cb94ee6SHong Zhang Level: advanced 6207cb94ee6SHong Zhang 6217cb94ee6SHong Zhang .keywords: Sundials, orthogonalization 6227cb94ee6SHong Zhang 6237cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6247cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 6257cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6267cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 6277cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 6287cb94ee6SHong Zhang 6297cb94ee6SHong Zhang @*/ 6307cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 6317cb94ee6SHong Zhang { 6327cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,TSSundialsGramSchmidtType); 6337cb94ee6SHong Zhang 6347cb94ee6SHong Zhang PetscFunctionBegin; 6357cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",(void (**)(void))&f);CHKERRQ(ierr); 6367cb94ee6SHong Zhang if (f) { 6377cb94ee6SHong Zhang ierr = (*f)(ts,type);CHKERRQ(ierr); 6387cb94ee6SHong Zhang } 6397cb94ee6SHong Zhang PetscFunctionReturn(0); 6407cb94ee6SHong Zhang } 6417cb94ee6SHong Zhang 6427cb94ee6SHong Zhang #undef __FUNCT__ 6437cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance" 6447cb94ee6SHong Zhang /*@ 6457cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 6467cb94ee6SHong Zhang Sundials for error control. 6477cb94ee6SHong Zhang 6487cb94ee6SHong Zhang Collective on TS 6497cb94ee6SHong Zhang 6507cb94ee6SHong Zhang Input parameters: 6517cb94ee6SHong Zhang + ts - the time-step context 6527cb94ee6SHong Zhang . aabs - the absolute tolerance 6537cb94ee6SHong Zhang - rel - the relative tolerance 6547cb94ee6SHong Zhang 6557cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 6567cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 6577cb94ee6SHong Zhang 6587cb94ee6SHong Zhang Level: intermediate 6597cb94ee6SHong Zhang 6607cb94ee6SHong Zhang .keywords: Sundials, tolerance 6617cb94ee6SHong Zhang 6627cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6637cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 6647cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6657cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 6667cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 6677cb94ee6SHong Zhang 6687cb94ee6SHong Zhang @*/ 6697cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance(TS ts,double aabs,double rel) 6707cb94ee6SHong Zhang { 6717cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,double,double); 6727cb94ee6SHong Zhang 6737cb94ee6SHong Zhang PetscFunctionBegin; 6747cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 6757cb94ee6SHong Zhang if (f) { 6767cb94ee6SHong Zhang ierr = (*f)(ts,aabs,rel);CHKERRQ(ierr); 6777cb94ee6SHong Zhang } 6787cb94ee6SHong Zhang PetscFunctionReturn(0); 6797cb94ee6SHong Zhang } 6807cb94ee6SHong Zhang 6817cb94ee6SHong Zhang #undef __FUNCT__ 6827cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC" 6837cb94ee6SHong Zhang /*@ 6847cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 6857cb94ee6SHong Zhang 6867cb94ee6SHong Zhang Input Parameter: 6877cb94ee6SHong Zhang . ts - the time-step context 6887cb94ee6SHong Zhang 6897cb94ee6SHong Zhang Output Parameter: 6907cb94ee6SHong Zhang . pc - the preconditioner context 6917cb94ee6SHong Zhang 6927cb94ee6SHong Zhang Level: advanced 6937cb94ee6SHong Zhang 6947cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6957cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 6967cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6977cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 6987cb94ee6SHong Zhang @*/ 6997cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC(TS ts,PC *pc) 7007cb94ee6SHong Zhang { 7017cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,PC *); 7027cb94ee6SHong Zhang 7037cb94ee6SHong Zhang PetscFunctionBegin; 7047cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetPC_C",(void (**)(void))&f);CHKERRQ(ierr); 7057cb94ee6SHong Zhang if (f) { 7067cb94ee6SHong Zhang ierr = (*f)(ts,pc);CHKERRQ(ierr); 7077cb94ee6SHong Zhang } else { 7087cb94ee6SHong Zhang SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"TS must be of Sundials type to extract the PC"); 7097cb94ee6SHong Zhang } 7107cb94ee6SHong Zhang 7117cb94ee6SHong Zhang PetscFunctionReturn(0); 7127cb94ee6SHong Zhang } 7137cb94ee6SHong Zhang 7147cb94ee6SHong Zhang #undef __FUNCT__ 7157cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime" 7167cb94ee6SHong Zhang /*@ 7177cb94ee6SHong Zhang TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the 7187cb94ee6SHong Zhang exact final time requested by the user or just returns it at the final time 7197cb94ee6SHong Zhang it computed. (Defaults to true). 7207cb94ee6SHong Zhang 7217cb94ee6SHong Zhang Input Parameter: 7227cb94ee6SHong Zhang + ts - the time-step context 7237cb94ee6SHong Zhang - ft - PETSC_TRUE if interpolates, else PETSC_FALSE 7247cb94ee6SHong Zhang 7257cb94ee6SHong Zhang Level: beginner 7267cb94ee6SHong Zhang 7277cb94ee6SHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7287cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 7297cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7307cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 7317cb94ee6SHong Zhang @*/ 7327cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime(TS ts,PetscTruth ft) 7337cb94ee6SHong Zhang { 7347cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,PetscTruth); 7357cb94ee6SHong Zhang 7367cb94ee6SHong Zhang PetscFunctionBegin; 7377cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetExactFinalTime_C",(void (**)(void))&f);CHKERRQ(ierr); 7387cb94ee6SHong Zhang if (f) { 7397cb94ee6SHong Zhang ierr = (*f)(ts,ft);CHKERRQ(ierr); 7407cb94ee6SHong Zhang } 7417cb94ee6SHong Zhang PetscFunctionReturn(0); 7427cb94ee6SHong Zhang } 7437cb94ee6SHong Zhang 7447cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 7457cb94ee6SHong Zhang /*MC 7467cb94ee6SHong Zhang TS_Sundials - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 7477cb94ee6SHong Zhang 7487cb94ee6SHong Zhang Options Database: 7497cb94ee6SHong Zhang + -ts_sundials_type <bdf,adams> 7507cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 7517cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 7527cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 7537cb94ee6SHong Zhang . -ts_sundials_linear_tolerance <tol> 7547cb94ee6SHong Zhang . -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions 7557cb94ee6SHong Zhang - -ts_sundials_not_exact_final_time -Allow SUNDIALS to stop near the final time, not exactly on it 7567cb94ee6SHong Zhang 7577cb94ee6SHong Zhang Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply 7587cb94ee6SHong Zhang only PETSc PC options 7597cb94ee6SHong Zhang 7607cb94ee6SHong Zhang Level: beginner 7617cb94ee6SHong Zhang 7627cb94ee6SHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(), 7637cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime() 7647cb94ee6SHong Zhang 7657cb94ee6SHong Zhang M*/ 7667cb94ee6SHong Zhang EXTERN_C_BEGIN 7677cb94ee6SHong Zhang #undef __FUNCT__ 7687cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials" 7697cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Sundials(TS ts) 7707cb94ee6SHong Zhang { 7717cb94ee6SHong Zhang TS_Sundials *cvode; 7727cb94ee6SHong Zhang PetscErrorCode ierr; 7737cb94ee6SHong Zhang 7747cb94ee6SHong Zhang PetscFunctionBegin; 7757cb94ee6SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 7767cb94ee6SHong Zhang ts->ops->view = TSView_Sundials; 7777cb94ee6SHong Zhang 7787cb94ee6SHong Zhang if (ts->problem_type != TS_NONLINEAR) { 7797cb94ee6SHong Zhang SETERRQ(PETSC_ERR_SUP,"Only support for nonlinear problems"); 7807cb94ee6SHong Zhang } 7817cb94ee6SHong Zhang ts->ops->setup = TSSetUp_Sundials_Nonlinear; 7827cb94ee6SHong Zhang ts->ops->step = TSStep_Sundials_Nonlinear; 7837cb94ee6SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_Sundials_Nonlinear; 7847cb94ee6SHong Zhang 7857cb94ee6SHong Zhang ierr = PetscNew(TS_Sundials,&cvode);CHKERRQ(ierr); 7867cb94ee6SHong Zhang ierr = PCCreate(ts->comm,&cvode->pc);CHKERRQ(ierr); 7877cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->pc);CHKERRQ(ierr); 7887cb94ee6SHong Zhang ts->data = (void*)cvode; 7897cb94ee6SHong Zhang cvode->cvode_type = CV_BDF; 7907cb94ee6SHong Zhang cvode->gtype = SUNDIALS_UNMODIFIED_GS; 7917cb94ee6SHong Zhang cvode->restart = 5; 7927cb94ee6SHong Zhang cvode->linear_tol = .05; 7937cb94ee6SHong Zhang 7947cb94ee6SHong Zhang cvode->exact_final_time = PETSC_FALSE; 7957cb94ee6SHong Zhang 7967cb94ee6SHong Zhang ierr = MPI_Comm_dup(ts->comm,&(cvode->comm_sundials));CHKERRQ(ierr); 7977cb94ee6SHong Zhang /* set tolerance for Sundials */ 7987cb94ee6SHong Zhang cvode->abstol = 1e-6; 7997cb94ee6SHong Zhang cvode->reltol = 1e-6; 8007cb94ee6SHong Zhang 8017cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials", 8027cb94ee6SHong Zhang TSSundialsSetType_Sundials);CHKERRQ(ierr); 8037cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C", 8047cb94ee6SHong Zhang "TSSundialsSetGMRESRestart_Sundials", 8057cb94ee6SHong Zhang TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr); 8067cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C", 8077cb94ee6SHong Zhang "TSSundialsSetLinearTolerance_Sundials", 8087cb94ee6SHong Zhang TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 8097cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C", 8107cb94ee6SHong Zhang "TSSundialsSetGramSchmidtType_Sundials", 8117cb94ee6SHong Zhang TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 8127cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C", 8137cb94ee6SHong Zhang "TSSundialsSetTolerance_Sundials", 8147cb94ee6SHong Zhang TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 8157cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C", 8167cb94ee6SHong Zhang "TSSundialsGetPC_Sundials", 8177cb94ee6SHong Zhang TSSundialsGetPC_Sundials);CHKERRQ(ierr); 8187cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C", 8197cb94ee6SHong Zhang "TSSundialsGetIterations_Sundials", 8207cb94ee6SHong Zhang TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 8217cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C", 8227cb94ee6SHong Zhang "TSSundialsSetExactFinalTime_Sundials", 8237cb94ee6SHong Zhang TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr); 8247cb94ee6SHong Zhang PetscFunctionReturn(0); 8257cb94ee6SHong Zhang } 8267cb94ee6SHong Zhang EXTERN_C_END 8277cb94ee6SHong Zhang 8287cb94ee6SHong Zhang 8297cb94ee6SHong Zhang 8307cb94ee6SHong Zhang 8317cb94ee6SHong Zhang 8327cb94ee6SHong Zhang 8337cb94ee6SHong Zhang 8347cb94ee6SHong Zhang 8357cb94ee6SHong Zhang 8367cb94ee6SHong Zhang 837