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 PetscErrorCode TSStep_Sundials_Nonlinear(TS ts,int *steps,double *time) 1327cb94ee6SHong Zhang { 1337cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1347cb94ee6SHong Zhang Vec sol = ts->vec_sol; 1357cb94ee6SHong Zhang PetscErrorCode ierr; 1365d47aee6SHong Zhang PetscInt i,max_steps = ts->max_steps,flag; 1377cb94ee6SHong Zhang long int its; 1387cb94ee6SHong Zhang realtype t,tout; 1397cb94ee6SHong Zhang PetscScalar *y_data; 1407cb94ee6SHong Zhang void *mem; 1417cb94ee6SHong Zhang 1427cb94ee6SHong Zhang PetscFunctionBegin; 14316016685SHong Zhang mem = cvode->mem; 1447cb94ee6SHong Zhang tout = ts->max_time; 1457cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr); 1467cb94ee6SHong Zhang N_VSetArrayPointer((realtype *)y_data,cvode->y); 1477cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 1487cb94ee6SHong Zhang for (i = 0; i < max_steps; i++) { 1497cb94ee6SHong Zhang if (ts->ptime >= ts->max_time) break; 1503f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 1512bfc04deSHong Zhang if (cvode->monitorstep){ 15228adb3f7SHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP); 1532bfc04deSHong Zhang } else { 1542bfc04deSHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL); 1552bfc04deSHong Zhang } 156e32f2f54SBarry Smith if (flag)SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag); 1577cb94ee6SHong Zhang if (t > ts->max_time && cvode->exact_final_time) { 1587cb94ee6SHong Zhang /* interpolate to final requested time */ 1597cb94ee6SHong Zhang ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr); 1607cb94ee6SHong Zhang t = tout; 1617cb94ee6SHong Zhang } 1627cb94ee6SHong Zhang ts->time_step = t - ts->ptime; 1637cb94ee6SHong Zhang ts->ptime = t; 1647cb94ee6SHong Zhang 1657cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 1667cb94ee6SHong Zhang ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr); 1677cb94ee6SHong Zhang ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr); 1687cb94ee6SHong Zhang ierr = VecResetArray(cvode->w1); CHKERRQ(ierr); 1697cb94ee6SHong Zhang ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr); 1707cb94ee6SHong Zhang ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr); 1717cb94ee6SHong Zhang ts->nonlinear_its = its; 1724ac7836bSHong Zhang ierr = CVSpilsGetNumLinIters(mem, &its); 1737cb94ee6SHong Zhang ts->linear_its = its; 1747cb94ee6SHong Zhang ts->steps++; 1753f2090d5SJed Brown ierr = TSPostStep(ts);CHKERRQ(ierr); 1767cb94ee6SHong Zhang ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr); 1777cb94ee6SHong Zhang } 1787cb94ee6SHong Zhang *steps += ts->steps; 1797cb94ee6SHong Zhang *time = t; 1807cb94ee6SHong Zhang PetscFunctionReturn(0); 1817cb94ee6SHong Zhang } 1827cb94ee6SHong Zhang 1837cb94ee6SHong Zhang #undef __FUNCT__ 1847cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials" 1857cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts) 1867cb94ee6SHong Zhang { 1877cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1887cb94ee6SHong Zhang PetscErrorCode ierr; 1897cb94ee6SHong Zhang 1907cb94ee6SHong Zhang PetscFunctionBegin; 1917cb94ee6SHong Zhang if (cvode->pmat) {ierr = MatDestroy(cvode->pmat);CHKERRQ(ierr);} 1927cb94ee6SHong Zhang if (cvode->pc) {ierr = PCDestroy(cvode->pc);CHKERRQ(ierr);} 1937cb94ee6SHong Zhang if (cvode->update) {ierr = VecDestroy(cvode->update);CHKERRQ(ierr);} 1947cb94ee6SHong Zhang if (cvode->func) {ierr = VecDestroy(cvode->func);CHKERRQ(ierr);} 1957cb94ee6SHong Zhang if (cvode->rhs) {ierr = VecDestroy(cvode->rhs);CHKERRQ(ierr);} 1967cb94ee6SHong Zhang if (cvode->w1) {ierr = VecDestroy(cvode->w1);CHKERRQ(ierr);} 1977cb94ee6SHong Zhang if (cvode->w2) {ierr = VecDestroy(cvode->w2);CHKERRQ(ierr);} 198a07356b8SSatish Balay ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr); 1992c823083SHong Zhang if (cvode->mem) {CVodeFree(&cvode->mem);} 2007cb94ee6SHong Zhang ierr = PetscFree(cvode);CHKERRQ(ierr); 2017cb94ee6SHong Zhang PetscFunctionReturn(0); 2027cb94ee6SHong Zhang } 2037cb94ee6SHong Zhang 2047cb94ee6SHong Zhang #undef __FUNCT__ 2057cb94ee6SHong Zhang #define __FUNCT__ "TSSetUp_Sundials_Nonlinear" 2067cb94ee6SHong Zhang PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts) 2077cb94ee6SHong Zhang { 2087cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2097cb94ee6SHong Zhang PetscErrorCode ierr; 21016016685SHong Zhang PetscInt glosize,locsize,i,flag; 2117cb94ee6SHong Zhang PetscScalar *y_data,*parray; 21216016685SHong Zhang void *mem; 21316016685SHong Zhang const PCType pctype; 21416016685SHong Zhang PetscTruth pcnone; 21516016685SHong Zhang Vec sol = ts->vec_sol; 2167cb94ee6SHong Zhang 2177cb94ee6SHong Zhang PetscFunctionBegin; 2187cb94ee6SHong Zhang ierr = PCSetFromOptions(cvode->pc);CHKERRQ(ierr); 2197cb94ee6SHong Zhang /* get the vector size */ 2207cb94ee6SHong Zhang ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr); 2217cb94ee6SHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr); 2227cb94ee6SHong Zhang 2237cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 224a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 225e32f2f54SBarry Smith if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated"); 2267cb94ee6SHong Zhang 22728adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 2287cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr); 2297cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y); 2307cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 2317cb94ee6SHong Zhang /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/ 2327cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 2337cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr); 2347cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->func);CHKERRQ(ierr); 2357cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr); 2367cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->func);CHKERRQ(ierr); 2377cb94ee6SHong Zhang 2387cb94ee6SHong Zhang /* 2397cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 2407cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 2417cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 2427cb94ee6SHong Zhang */ 2437adad957SLisandro Dalcin ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr); 2447adad957SLisandro Dalcin ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr); 2457cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr); 2467cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr); 24716016685SHong Zhang 24816016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 24916016685SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 250e32f2f54SBarry Smith if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails"); 25116016685SHong Zhang cvode->mem = mem; 25216016685SHong Zhang 25316016685SHong Zhang /* Set the pointer to user-defined data */ 25416016685SHong Zhang flag = CVodeSetUserData(mem, ts); 255e32f2f54SBarry Smith if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails"); 25616016685SHong Zhang 257*fc6b9e64SJed Brown /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 258*fc6b9e64SJed Brown flag = CVodeSetInitStep(mem,(realtype)ts->initial_time_step); 259*fc6b9e64SJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed"); 260f1cd61daSJed Brown if (cvode->mindt > 0) { 261f1cd61daSJed Brown flag = CVodeSetMinStep(mem,(realtype)cvode->mindt); 262f1cd61daSJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed"); 263f1cd61daSJed Brown } 264f1cd61daSJed Brown if (cvode->maxdt > 0) { 265f1cd61daSJed Brown flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt); 266f1cd61daSJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed"); 267f1cd61daSJed Brown } 268f1cd61daSJed Brown 26916016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 27016016685SHong Zhang * user's right hand side function in u'=f(t,u), the inital time T0, and 27116016685SHong Zhang * the initial dependent variable vector cvode->y */ 27216016685SHong Zhang flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 27316016685SHong Zhang if (flag){ 274e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag); 27516016685SHong Zhang } 27616016685SHong Zhang 27716016685SHong Zhang flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 27816016685SHong Zhang if (flag){ 279e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag); 28016016685SHong Zhang } 28116016685SHong Zhang 28216016685SHong Zhang /* initialize the number of steps */ 28316016685SHong Zhang ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 28416016685SHong Zhang 28516016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 28616016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 28716016685SHong Zhang ierr = PCGetType(cvode->pc,&pctype);CHKERRQ(ierr); 28816016685SHong Zhang ierr = PetscTypeCompare((PetscObject)cvode->pc,PCNONE,&pcnone);CHKERRQ(ierr); 28916016685SHong Zhang if (pcnone){ 29016016685SHong Zhang flag = CVSpgmr(mem,PREC_NONE,0); 291e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 29216016685SHong Zhang } else { 29316016685SHong Zhang flag = CVSpgmr(mem,PREC_LEFT,0); 294e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 29516016685SHong Zhang 29616016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 29716016685SHong Zhang and the pointer to the user-defined block data */ 29816016685SHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 299e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 30016016685SHong Zhang } 30116016685SHong Zhang 30216016685SHong Zhang flag = CVSpilsSetGSType(mem, MODIFIED_GS); 303e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag); 3047cb94ee6SHong Zhang PetscFunctionReturn(0); 3057cb94ee6SHong Zhang } 3067cb94ee6SHong Zhang 3076fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 308dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0}; 3096fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 310dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0}; 311a04cf4d8SBarry Smith 3127cb94ee6SHong Zhang #undef __FUNCT__ 3137cb94ee6SHong Zhang #define __FUNCT__ "TSSetFromOptions_Sundials_Nonlinear" 3147cb94ee6SHong Zhang PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts) 3157cb94ee6SHong Zhang { 3167cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 3177cb94ee6SHong Zhang PetscErrorCode ierr; 3187cb94ee6SHong Zhang int indx; 3197cb94ee6SHong Zhang PetscTruth flag; 3207cb94ee6SHong Zhang 3217cb94ee6SHong Zhang PetscFunctionBegin; 3227cb94ee6SHong Zhang ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr); 3236fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr); 3247cb94ee6SHong Zhang if (flag) { 3256fadb2cdSHong Zhang ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr); 3267cb94ee6SHong Zhang } 3276fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr); 3287cb94ee6SHong Zhang if (flag) { 3297cb94ee6SHong Zhang ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 3307cb94ee6SHong Zhang } 3317cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr); 3327cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr); 333f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr); 334f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr); 3357cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr); 3367cb94ee6SHong Zhang ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr); 33716016685SHong Zhang ierr = PetscOptionsTruth("-ts_sundials_exact_final_time","Interpolate output to stop exactly at the final time","TSSundialsSetExactFinalTime",cvode->exact_final_time,&cvode->exact_final_time,PETSC_NULL);CHKERRQ(ierr); 3382bfc04deSHong Zhang ierr = PetscOptionsTruth("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr); 3397cb94ee6SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 3407cb94ee6SHong Zhang PetscFunctionReturn(0); 3417cb94ee6SHong Zhang } 3427cb94ee6SHong Zhang 3437cb94ee6SHong Zhang #undef __FUNCT__ 3447cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials" 3457cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 3467cb94ee6SHong Zhang { 3477cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 3487cb94ee6SHong Zhang PetscErrorCode ierr; 3497cb94ee6SHong Zhang char *type; 3507cb94ee6SHong Zhang char atype[] = "Adams"; 3517cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 3527cb94ee6SHong Zhang PetscTruth iascii,isstring; 3532c823083SHong Zhang long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 3542c823083SHong Zhang PetscInt qlast,qcur; 3555d47aee6SHong Zhang PetscReal hinused,hlast,hcur,tcur,tolsfac; 3567cb94ee6SHong Zhang 3577cb94ee6SHong Zhang PetscFunctionBegin; 3587cb94ee6SHong Zhang if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;} 3597cb94ee6SHong Zhang else {type = btype;} 3607cb94ee6SHong Zhang 3612692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 3622692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 3637cb94ee6SHong Zhang if (iascii) { 3647cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr); 3657cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr); 3667cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr); 3677cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr); 3687cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr); 3697cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 3707cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 3717cb94ee6SHong Zhang } else { 3727cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 3737cb94ee6SHong Zhang } 374f1cd61daSJed Brown if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);} 375f1cd61daSJed Brown if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);} 3762c823083SHong Zhang 3775d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 3785d47aee6SHong Zhang ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr); 3795d47aee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr); 3802c823083SHong Zhang ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals, 3812c823083SHong Zhang &nlinsetups,&nfails,&qlast,&qcur, 3822c823083SHong Zhang &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr); 3832c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr); 3842c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr); 3852c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr); 3862c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr); 3872c823083SHong Zhang 3882c823083SHong Zhang ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr); 3892c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr); 3902c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr); 3912c823083SHong Zhang 3922c823083SHong Zhang ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */ 3932c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr); 3942c823083SHong Zhang ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr); 3952c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr); 3962c823083SHong Zhang ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr); 3972c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr); 3982c823083SHong Zhang ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr); 3992c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr); 4002c823083SHong Zhang ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4012c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr); 4022c823083SHong Zhang ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4032c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr); 4047cb94ee6SHong Zhang } else if (isstring) { 4057cb94ee6SHong Zhang ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 4067cb94ee6SHong Zhang } else { 407e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name); 4087cb94ee6SHong Zhang } 4097cb94ee6SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 4107cb94ee6SHong Zhang ierr = PCView(cvode->pc,viewer);CHKERRQ(ierr); 4117cb94ee6SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 4127cb94ee6SHong Zhang PetscFunctionReturn(0); 4137cb94ee6SHong Zhang } 4147cb94ee6SHong Zhang 4157cb94ee6SHong Zhang 4167cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 4177cb94ee6SHong Zhang EXTERN_C_BEGIN 4187cb94ee6SHong Zhang #undef __FUNCT__ 4197cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials" 4206fadb2cdSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 4217cb94ee6SHong Zhang { 4227cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4237cb94ee6SHong Zhang 4247cb94ee6SHong Zhang PetscFunctionBegin; 4257cb94ee6SHong Zhang cvode->cvode_type = type; 4267cb94ee6SHong Zhang PetscFunctionReturn(0); 4277cb94ee6SHong Zhang } 4287cb94ee6SHong Zhang EXTERN_C_END 4297cb94ee6SHong Zhang 4307cb94ee6SHong Zhang EXTERN_C_BEGIN 4317cb94ee6SHong Zhang #undef __FUNCT__ 4327cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials" 4337cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart_Sundials(TS ts,int restart) 4347cb94ee6SHong Zhang { 4357cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4367cb94ee6SHong Zhang 4377cb94ee6SHong Zhang PetscFunctionBegin; 4387cb94ee6SHong Zhang cvode->restart = restart; 4397cb94ee6SHong Zhang PetscFunctionReturn(0); 4407cb94ee6SHong Zhang } 4417cb94ee6SHong Zhang EXTERN_C_END 4427cb94ee6SHong Zhang 4437cb94ee6SHong Zhang EXTERN_C_BEGIN 4447cb94ee6SHong Zhang #undef __FUNCT__ 4457cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials" 4467cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance_Sundials(TS ts,double tol) 4477cb94ee6SHong Zhang { 4487cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4497cb94ee6SHong Zhang 4507cb94ee6SHong Zhang PetscFunctionBegin; 4517cb94ee6SHong Zhang cvode->linear_tol = tol; 4527cb94ee6SHong Zhang PetscFunctionReturn(0); 4537cb94ee6SHong Zhang } 4547cb94ee6SHong Zhang EXTERN_C_END 4557cb94ee6SHong Zhang 4567cb94ee6SHong Zhang EXTERN_C_BEGIN 4577cb94ee6SHong Zhang #undef __FUNCT__ 4587cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials" 4597cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 4607cb94ee6SHong Zhang { 4617cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4627cb94ee6SHong Zhang 4637cb94ee6SHong Zhang PetscFunctionBegin; 4647cb94ee6SHong Zhang cvode->gtype = type; 4657cb94ee6SHong Zhang PetscFunctionReturn(0); 4667cb94ee6SHong Zhang } 4677cb94ee6SHong Zhang EXTERN_C_END 4687cb94ee6SHong Zhang 4697cb94ee6SHong Zhang EXTERN_C_BEGIN 4707cb94ee6SHong Zhang #undef __FUNCT__ 4717cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials" 4727cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel) 4737cb94ee6SHong Zhang { 4747cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4757cb94ee6SHong Zhang 4767cb94ee6SHong Zhang PetscFunctionBegin; 4777cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 4787cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 4797cb94ee6SHong Zhang PetscFunctionReturn(0); 4807cb94ee6SHong Zhang } 4817cb94ee6SHong Zhang EXTERN_C_END 4827cb94ee6SHong Zhang 4837cb94ee6SHong Zhang EXTERN_C_BEGIN 4847cb94ee6SHong Zhang #undef __FUNCT__ 485f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials" 486f1cd61daSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 487f1cd61daSJed Brown { 488f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 489f1cd61daSJed Brown 490f1cd61daSJed Brown PetscFunctionBegin; 491f1cd61daSJed Brown cvode->mindt = mindt; 492f1cd61daSJed Brown PetscFunctionReturn(0); 493f1cd61daSJed Brown } 494f1cd61daSJed Brown EXTERN_C_END 495f1cd61daSJed Brown 496f1cd61daSJed Brown EXTERN_C_BEGIN 497f1cd61daSJed Brown #undef __FUNCT__ 498f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials" 499f1cd61daSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 500f1cd61daSJed Brown { 501f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 502f1cd61daSJed Brown 503f1cd61daSJed Brown PetscFunctionBegin; 504f1cd61daSJed Brown cvode->maxdt = maxdt; 505f1cd61daSJed Brown PetscFunctionReturn(0); 506f1cd61daSJed Brown } 507f1cd61daSJed Brown EXTERN_C_END 508f1cd61daSJed Brown 509f1cd61daSJed Brown EXTERN_C_BEGIN 510f1cd61daSJed Brown #undef __FUNCT__ 5117cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials" 5127cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC_Sundials(TS ts,PC *pc) 5137cb94ee6SHong Zhang { 5147cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5157cb94ee6SHong Zhang 5167cb94ee6SHong Zhang PetscFunctionBegin; 5177cb94ee6SHong Zhang *pc = cvode->pc; 5187cb94ee6SHong Zhang PetscFunctionReturn(0); 5197cb94ee6SHong Zhang } 5207cb94ee6SHong Zhang EXTERN_C_END 5217cb94ee6SHong Zhang 5227cb94ee6SHong Zhang EXTERN_C_BEGIN 5237cb94ee6SHong Zhang #undef __FUNCT__ 5247cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials" 5257cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 5267cb94ee6SHong Zhang { 5277cb94ee6SHong Zhang PetscFunctionBegin; 5282c823083SHong Zhang if (nonlin) *nonlin = ts->nonlinear_its; 5292c823083SHong Zhang if (lin) *lin = ts->linear_its; 5307cb94ee6SHong Zhang PetscFunctionReturn(0); 5317cb94ee6SHong Zhang } 5327cb94ee6SHong Zhang EXTERN_C_END 5337cb94ee6SHong Zhang 5347cb94ee6SHong Zhang EXTERN_C_BEGIN 5357cb94ee6SHong Zhang #undef __FUNCT__ 5367cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials" 5377cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime_Sundials(TS ts,PetscTruth s) 5387cb94ee6SHong Zhang { 5397cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5407cb94ee6SHong Zhang 5417cb94ee6SHong Zhang PetscFunctionBegin; 5427cb94ee6SHong Zhang cvode->exact_final_time = s; 5437cb94ee6SHong Zhang PetscFunctionReturn(0); 5447cb94ee6SHong Zhang } 5457cb94ee6SHong Zhang EXTERN_C_END 5462bfc04deSHong Zhang 5472bfc04deSHong Zhang EXTERN_C_BEGIN 5482bfc04deSHong Zhang #undef __FUNCT__ 5492bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials" 5502bfc04deSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscTruth s) 5512bfc04deSHong Zhang { 5522bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5532bfc04deSHong Zhang 5542bfc04deSHong Zhang PetscFunctionBegin; 5552bfc04deSHong Zhang cvode->monitorstep = s; 5562bfc04deSHong Zhang PetscFunctionReturn(0); 5572bfc04deSHong Zhang } 5582bfc04deSHong Zhang EXTERN_C_END 5597cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 5607cb94ee6SHong Zhang 5617cb94ee6SHong Zhang #undef __FUNCT__ 5627cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations" 5637cb94ee6SHong Zhang /*@C 5647cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 5657cb94ee6SHong Zhang 5667cb94ee6SHong Zhang Not Collective 5677cb94ee6SHong Zhang 5687cb94ee6SHong Zhang Input parameters: 5697cb94ee6SHong Zhang . ts - the time-step context 5707cb94ee6SHong Zhang 5717cb94ee6SHong Zhang Output Parameters: 5727cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 5737cb94ee6SHong Zhang - lin - number of linear iterations 5747cb94ee6SHong Zhang 5757cb94ee6SHong Zhang Level: advanced 5767cb94ee6SHong Zhang 5777cb94ee6SHong Zhang Notes: 5787cb94ee6SHong Zhang These return the number since the creation of the TS object 5797cb94ee6SHong Zhang 5807cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations 5817cb94ee6SHong Zhang 5827cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(), 5837cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 5847cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 585f1cd61daSJed Brown TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSundialsSetExactFinalTime() 5867cb94ee6SHong Zhang 5877cb94ee6SHong Zhang @*/ 5887cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 5897cb94ee6SHong Zhang { 5907cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,int*,int*); 5917cb94ee6SHong Zhang 5927cb94ee6SHong Zhang PetscFunctionBegin; 5937cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetIterations_C",(void (**)(void))&f);CHKERRQ(ierr); 5947cb94ee6SHong Zhang if (f) { 5957cb94ee6SHong Zhang ierr = (*f)(ts,nonlin,lin);CHKERRQ(ierr); 5967cb94ee6SHong Zhang } 5977cb94ee6SHong Zhang PetscFunctionReturn(0); 5987cb94ee6SHong Zhang } 5997cb94ee6SHong Zhang 6007cb94ee6SHong Zhang #undef __FUNCT__ 6017cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType" 6027cb94ee6SHong Zhang /*@ 6037cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 6047cb94ee6SHong Zhang 6053f9fe445SBarry Smith Logically Collective on TS 6067cb94ee6SHong Zhang 6077cb94ee6SHong Zhang Input parameters: 6087cb94ee6SHong Zhang + ts - the time-step context 6097cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 6107cb94ee6SHong Zhang 6117cb94ee6SHong Zhang Level: intermediate 6127cb94ee6SHong Zhang 6137cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula 6147cb94ee6SHong Zhang 6157cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetGMRESRestart(), 6167cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 6177cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6187cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 6197cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 6207cb94ee6SHong Zhang @*/ 6216fadb2cdSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType(TS ts,TSSundialsLmmType type) 6227cb94ee6SHong Zhang { 6236fadb2cdSHong Zhang PetscErrorCode ierr,(*f)(TS,TSSundialsLmmType); 6247cb94ee6SHong Zhang 6257cb94ee6SHong Zhang PetscFunctionBegin; 6267cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 6277cb94ee6SHong Zhang if (f) { 6287cb94ee6SHong Zhang ierr = (*f)(ts,type);CHKERRQ(ierr); 6297cb94ee6SHong Zhang } 6307cb94ee6SHong Zhang PetscFunctionReturn(0); 6317cb94ee6SHong Zhang } 6327cb94ee6SHong Zhang 6337cb94ee6SHong Zhang #undef __FUNCT__ 6347cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart" 6357cb94ee6SHong Zhang /*@ 6367cb94ee6SHong Zhang TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by 6377cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 6387cb94ee6SHong Zhang this is ALSO the maximum number of GMRES steps that will be used. 6397cb94ee6SHong Zhang 6403f9fe445SBarry Smith Logically Collective on TS 6417cb94ee6SHong Zhang 6427cb94ee6SHong Zhang Input parameters: 6437cb94ee6SHong Zhang + ts - the time-step context 6447cb94ee6SHong Zhang - restart - number of direction vectors (the restart size). 6457cb94ee6SHong Zhang 6467cb94ee6SHong Zhang Level: advanced 6477cb94ee6SHong Zhang 6487cb94ee6SHong Zhang .keywords: GMRES, restart 6497cb94ee6SHong Zhang 6507cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 6517cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 6527cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6537cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 6547cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 6557cb94ee6SHong Zhang 6567cb94ee6SHong Zhang @*/ 6577cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart(TS ts,int restart) 6587cb94ee6SHong Zhang { 6597cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,int); 6607cb94ee6SHong Zhang 6617cb94ee6SHong Zhang PetscFunctionBegin; 662c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(ts,restart,2); 6637cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGMRESRestart_C",(void (**)(void))&f);CHKERRQ(ierr); 6647cb94ee6SHong Zhang if (f) { 6657cb94ee6SHong Zhang ierr = (*f)(ts,restart);CHKERRQ(ierr); 6667cb94ee6SHong Zhang } 6677cb94ee6SHong Zhang PetscFunctionReturn(0); 6687cb94ee6SHong Zhang } 6697cb94ee6SHong Zhang 6707cb94ee6SHong Zhang #undef __FUNCT__ 6717cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance" 6727cb94ee6SHong Zhang /*@ 6737cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 6747cb94ee6SHong Zhang system by SUNDIALS. 6757cb94ee6SHong Zhang 6763f9fe445SBarry Smith Logically Collective on TS 6777cb94ee6SHong Zhang 6787cb94ee6SHong Zhang Input parameters: 6797cb94ee6SHong Zhang + ts - the time-step context 6807cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 6817cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 6827cb94ee6SHong Zhang 6837cb94ee6SHong Zhang Level: advanced 6847cb94ee6SHong Zhang 6857cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS 6867cb94ee6SHong Zhang 6877cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6887cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 6897cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6907cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 6917cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 6927cb94ee6SHong Zhang 6937cb94ee6SHong Zhang @*/ 6947cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance(TS ts,double tol) 6957cb94ee6SHong Zhang { 6967cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,double); 6977cb94ee6SHong Zhang 6987cb94ee6SHong Zhang PetscFunctionBegin; 699c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,tol,2); 7007cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 7017cb94ee6SHong Zhang if (f) { 7027cb94ee6SHong Zhang ierr = (*f)(ts,tol);CHKERRQ(ierr); 7037cb94ee6SHong Zhang } 7047cb94ee6SHong Zhang PetscFunctionReturn(0); 7057cb94ee6SHong Zhang } 7067cb94ee6SHong Zhang 7077cb94ee6SHong Zhang #undef __FUNCT__ 7087cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType" 7097cb94ee6SHong Zhang /*@ 7107cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 7117cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 7127cb94ee6SHong Zhang 7133f9fe445SBarry Smith Logically Collective on TS 7147cb94ee6SHong Zhang 7157cb94ee6SHong Zhang Input parameters: 7167cb94ee6SHong Zhang + ts - the time-step context 7177cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 7187cb94ee6SHong Zhang 7197cb94ee6SHong Zhang Level: advanced 7207cb94ee6SHong Zhang 7217cb94ee6SHong Zhang .keywords: Sundials, orthogonalization 7227cb94ee6SHong Zhang 7237cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7247cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 7257cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7267cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 7277cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 7287cb94ee6SHong Zhang 7297cb94ee6SHong Zhang @*/ 7307cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 7317cb94ee6SHong Zhang { 7327cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,TSSundialsGramSchmidtType); 7337cb94ee6SHong Zhang 7347cb94ee6SHong Zhang PetscFunctionBegin; 7357cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",(void (**)(void))&f);CHKERRQ(ierr); 7367cb94ee6SHong Zhang if (f) { 7377cb94ee6SHong Zhang ierr = (*f)(ts,type);CHKERRQ(ierr); 7387cb94ee6SHong Zhang } 7397cb94ee6SHong Zhang PetscFunctionReturn(0); 7407cb94ee6SHong Zhang } 7417cb94ee6SHong Zhang 7427cb94ee6SHong Zhang #undef __FUNCT__ 7437cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance" 7447cb94ee6SHong Zhang /*@ 7457cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 7467cb94ee6SHong Zhang Sundials for error control. 7477cb94ee6SHong Zhang 7483f9fe445SBarry Smith Logically Collective on TS 7497cb94ee6SHong Zhang 7507cb94ee6SHong Zhang Input parameters: 7517cb94ee6SHong Zhang + ts - the time-step context 7527cb94ee6SHong Zhang . aabs - the absolute tolerance 7537cb94ee6SHong Zhang - rel - the relative tolerance 7547cb94ee6SHong Zhang 7557cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 7567cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 7577cb94ee6SHong Zhang 7587cb94ee6SHong Zhang Level: intermediate 7597cb94ee6SHong Zhang 7607cb94ee6SHong Zhang .keywords: Sundials, tolerance 7617cb94ee6SHong Zhang 7627cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7637cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 7647cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7657cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 7667cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 7677cb94ee6SHong Zhang 7687cb94ee6SHong Zhang @*/ 7697cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance(TS ts,double aabs,double rel) 7707cb94ee6SHong Zhang { 7717cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,double,double); 7727cb94ee6SHong Zhang 7737cb94ee6SHong Zhang PetscFunctionBegin; 7747cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 7757cb94ee6SHong Zhang if (f) { 7767cb94ee6SHong Zhang ierr = (*f)(ts,aabs,rel);CHKERRQ(ierr); 7777cb94ee6SHong Zhang } 7787cb94ee6SHong Zhang PetscFunctionReturn(0); 7797cb94ee6SHong Zhang } 7807cb94ee6SHong Zhang 7817cb94ee6SHong Zhang #undef __FUNCT__ 7827cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC" 7837cb94ee6SHong Zhang /*@ 7847cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 7857cb94ee6SHong Zhang 7867cb94ee6SHong Zhang Input Parameter: 7877cb94ee6SHong Zhang . ts - the time-step context 7887cb94ee6SHong Zhang 7897cb94ee6SHong Zhang Output Parameter: 7907cb94ee6SHong Zhang . pc - the preconditioner context 7917cb94ee6SHong Zhang 7927cb94ee6SHong Zhang Level: advanced 7937cb94ee6SHong Zhang 7947cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7957cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 7967cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7977cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 7987cb94ee6SHong Zhang @*/ 7997cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC(TS ts,PC *pc) 8007cb94ee6SHong Zhang { 8017cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,PC *); 8027cb94ee6SHong Zhang 8037cb94ee6SHong Zhang PetscFunctionBegin; 8047cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetPC_C",(void (**)(void))&f);CHKERRQ(ierr); 8057cb94ee6SHong Zhang if (f) { 8067cb94ee6SHong Zhang ierr = (*f)(ts,pc);CHKERRQ(ierr); 8077cb94ee6SHong Zhang } else { 808e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TS must be of Sundials type to extract the PC"); 8097cb94ee6SHong Zhang } 8107cb94ee6SHong Zhang 8117cb94ee6SHong Zhang PetscFunctionReturn(0); 8127cb94ee6SHong Zhang } 8137cb94ee6SHong Zhang 8147cb94ee6SHong Zhang #undef __FUNCT__ 8157cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime" 8167cb94ee6SHong Zhang /*@ 8177cb94ee6SHong Zhang TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the 8187cb94ee6SHong Zhang exact final time requested by the user or just returns it at the final time 8197cb94ee6SHong Zhang it computed. (Defaults to true). 8207cb94ee6SHong Zhang 8217cb94ee6SHong Zhang Input Parameter: 8227cb94ee6SHong Zhang + ts - the time-step context 8237cb94ee6SHong Zhang - ft - PETSC_TRUE if interpolates, else PETSC_FALSE 8247cb94ee6SHong Zhang 8257cb94ee6SHong Zhang Level: beginner 8267cb94ee6SHong Zhang 8277cb94ee6SHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8287cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 8297cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8307cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 8317cb94ee6SHong Zhang @*/ 8327cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime(TS ts,PetscTruth ft) 8337cb94ee6SHong Zhang { 8347cb94ee6SHong Zhang PetscErrorCode ierr,(*f)(TS,PetscTruth); 8357cb94ee6SHong Zhang 8367cb94ee6SHong Zhang PetscFunctionBegin; 8377cb94ee6SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetExactFinalTime_C",(void (**)(void))&f);CHKERRQ(ierr); 8387cb94ee6SHong Zhang if (f) { 8397cb94ee6SHong Zhang ierr = (*f)(ts,ft);CHKERRQ(ierr); 8407cb94ee6SHong Zhang } 8417cb94ee6SHong Zhang PetscFunctionReturn(0); 8427cb94ee6SHong Zhang } 8437cb94ee6SHong Zhang 8442bfc04deSHong Zhang #undef __FUNCT__ 845f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep" 846f1cd61daSJed Brown /*@ 847f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 848f1cd61daSJed Brown 849f1cd61daSJed Brown Input Parameter: 850f1cd61daSJed Brown + ts - the time-step context 851f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 852f1cd61daSJed Brown 853*fc6b9e64SJed Brown Note: 854*fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 855*fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 856*fc6b9e64SJed Brown 857f1cd61daSJed Brown Level: beginner 858f1cd61daSJed Brown 859f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 860f1cd61daSJed Brown @*/ 861f1cd61daSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 862f1cd61daSJed Brown { 863f1cd61daSJed Brown PetscErrorCode ierr,(*f)(TS,PetscReal); 864f1cd61daSJed Brown 865f1cd61daSJed Brown PetscFunctionBegin; 866f1cd61daSJed Brown ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",(void(**)(void))&f);CHKERRQ(ierr); 867f1cd61daSJed Brown if (f) {ierr = (*f)(ts,mindt);CHKERRQ(ierr);} 868f1cd61daSJed Brown PetscFunctionReturn(0); 869f1cd61daSJed Brown } 870f1cd61daSJed Brown 871f1cd61daSJed Brown #undef __FUNCT__ 872f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep" 873f1cd61daSJed Brown /*@ 874f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 875f1cd61daSJed Brown 876f1cd61daSJed Brown Input Parameter: 877f1cd61daSJed Brown + ts - the time-step context 878f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 879f1cd61daSJed Brown 880f1cd61daSJed Brown Level: beginner 881f1cd61daSJed Brown 882f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 883f1cd61daSJed Brown @*/ 884f1cd61daSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 885f1cd61daSJed Brown { 886f1cd61daSJed Brown PetscErrorCode ierr,(*f)(TS,PetscReal); 887f1cd61daSJed Brown 888f1cd61daSJed Brown PetscFunctionBegin; 889f1cd61daSJed Brown ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",(void(**)(void))&f);CHKERRQ(ierr); 890f1cd61daSJed Brown if (f) {ierr = (*f)(ts,maxdt);CHKERRQ(ierr);} 891f1cd61daSJed Brown PetscFunctionReturn(0); 892f1cd61daSJed Brown } 893f1cd61daSJed Brown 894f1cd61daSJed Brown #undef __FUNCT__ 8952bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps" 8962bfc04deSHong Zhang /*@ 8972bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 8982bfc04deSHong Zhang 8992bfc04deSHong Zhang Input Parameter: 9002bfc04deSHong Zhang + ts - the time-step context 9012bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 9022bfc04deSHong Zhang 9032bfc04deSHong Zhang Level: beginner 9042bfc04deSHong Zhang 9052bfc04deSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 9062bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 9072bfc04deSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 9082bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 9092bfc04deSHong Zhang @*/ 9102bfc04deSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsMonitorInternalSteps(TS ts,PetscTruth ft) 9112bfc04deSHong Zhang { 9122bfc04deSHong Zhang PetscErrorCode ierr,(*f)(TS,PetscTruth); 9132bfc04deSHong Zhang 9142bfc04deSHong Zhang PetscFunctionBegin; 9152bfc04deSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",(void (**)(void))&f);CHKERRQ(ierr); 9162bfc04deSHong Zhang if (f) { 9172bfc04deSHong Zhang ierr = (*f)(ts,ft);CHKERRQ(ierr); 9182bfc04deSHong Zhang } 9192bfc04deSHong Zhang PetscFunctionReturn(0); 9202bfc04deSHong Zhang } 9217cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 9227cb94ee6SHong Zhang /*MC 92396f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 9247cb94ee6SHong Zhang 9257cb94ee6SHong Zhang Options Database: 9267cb94ee6SHong Zhang + -ts_sundials_type <bdf,adams> 9277cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 9287cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 9297cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 9307cb94ee6SHong Zhang . -ts_sundials_linear_tolerance <tol> 9317cb94ee6SHong Zhang . -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions 93216016685SHong Zhang . -ts_sundials_exact_final_time - Interpolate output to stop exactly at the final time 93316016685SHong Zhang - -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps 93416016685SHong Zhang 9357cb94ee6SHong Zhang 9367cb94ee6SHong Zhang Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply 9377cb94ee6SHong Zhang only PETSc PC options 9387cb94ee6SHong Zhang 9397cb94ee6SHong Zhang Level: beginner 9407cb94ee6SHong Zhang 9417cb94ee6SHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(), 9427cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime() 9437cb94ee6SHong Zhang 9447cb94ee6SHong Zhang M*/ 9457cb94ee6SHong Zhang EXTERN_C_BEGIN 9467cb94ee6SHong Zhang #undef __FUNCT__ 9477cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials" 9487cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Sundials(TS ts) 9497cb94ee6SHong Zhang { 9507cb94ee6SHong Zhang TS_Sundials *cvode; 9517cb94ee6SHong Zhang PetscErrorCode ierr; 9527cb94ee6SHong Zhang 9537cb94ee6SHong Zhang PetscFunctionBegin; 95417186662SBarry Smith if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for nonlinear problems"); 95528adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 95628adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 9577cb94ee6SHong Zhang ts->ops->setup = TSSetUp_Sundials_Nonlinear; 9587cb94ee6SHong Zhang ts->ops->step = TSStep_Sundials_Nonlinear; 9597cb94ee6SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_Sundials_Nonlinear; 9607cb94ee6SHong Zhang 96138f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr); 9627adad957SLisandro Dalcin ierr = PCCreate(((PetscObject)ts)->comm,&cvode->pc);CHKERRQ(ierr); 9637cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->pc);CHKERRQ(ierr); 9647cb94ee6SHong Zhang ts->data = (void*)cvode; 9656fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 9666fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 9677cb94ee6SHong Zhang cvode->restart = 5; 9687cb94ee6SHong Zhang cvode->linear_tol = .05; 9697cb94ee6SHong Zhang 9702bfc04deSHong Zhang cvode->exact_final_time = PETSC_TRUE; 9712bfc04deSHong Zhang cvode->monitorstep = PETSC_FALSE; 9727cb94ee6SHong Zhang 973a07356b8SSatish Balay ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr); 974f1cd61daSJed Brown 975f1cd61daSJed Brown cvode->mindt = -1.; 976f1cd61daSJed Brown cvode->maxdt = -1.; 977f1cd61daSJed Brown 9787cb94ee6SHong Zhang /* set tolerance for Sundials */ 9797cb94ee6SHong Zhang cvode->reltol = 1e-6; 9802c823083SHong Zhang cvode->abstol = 1e-6; 9817cb94ee6SHong Zhang 9827cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials", 9837cb94ee6SHong Zhang TSSundialsSetType_Sundials);CHKERRQ(ierr); 9847cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C", 9857cb94ee6SHong Zhang "TSSundialsSetGMRESRestart_Sundials", 9867cb94ee6SHong Zhang TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr); 9877cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C", 9887cb94ee6SHong Zhang "TSSundialsSetLinearTolerance_Sundials", 9897cb94ee6SHong Zhang TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 9907cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C", 9917cb94ee6SHong Zhang "TSSundialsSetGramSchmidtType_Sundials", 9927cb94ee6SHong Zhang TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 9937cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C", 9947cb94ee6SHong Zhang "TSSundialsSetTolerance_Sundials", 9957cb94ee6SHong Zhang TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 996f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C", 997f1cd61daSJed Brown "TSSundialsSetMinTimeStep_Sundials", 998f1cd61daSJed Brown TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr); 999f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C", 1000f1cd61daSJed Brown "TSSundialsSetMaxTimeStep_Sundials", 1001f1cd61daSJed Brown TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr); 10027cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C", 10037cb94ee6SHong Zhang "TSSundialsGetPC_Sundials", 10047cb94ee6SHong Zhang TSSundialsGetPC_Sundials);CHKERRQ(ierr); 10057cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C", 10067cb94ee6SHong Zhang "TSSundialsGetIterations_Sundials", 10077cb94ee6SHong Zhang TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 10087cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C", 10097cb94ee6SHong Zhang "TSSundialsSetExactFinalTime_Sundials", 10107cb94ee6SHong Zhang TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr); 10112bfc04deSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C", 10122bfc04deSHong Zhang "TSSundialsMonitorInternalSteps_Sundials", 10132bfc04deSHong Zhang TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr); 10142bfc04deSHong Zhang 10157cb94ee6SHong Zhang PetscFunctionReturn(0); 10167cb94ee6SHong Zhang } 10177cb94ee6SHong Zhang EXTERN_C_END 1018