17cb94ee6SHong Zhang 27cb94ee6SHong Zhang /* 3db268bfaSHong Zhang Provides a PETSc interface to SUNDIALS/CVODE solver. 47cb94ee6SHong Zhang The interface to PVODE (old version of CVODE) was originally contributed 57cb94ee6SHong Zhang by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik. 628adb3f7SHong Zhang 728adb3f7SHong Zhang Reference: sundials-2.4.0/examples/cvode/parallel/cvDiurnal_kry_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); 532c823083SHong Zhang } else { 547cb94ee6SHong Zhang ierr = MatCopy(Jac,cvode->pmat,str);CHKERRQ(ierr); 557cb94ee6SHong Zhang } 567cb94ee6SHong Zhang *jcurPtr = TRUE; 577cb94ee6SHong Zhang } 587cb94ee6SHong Zhang 597cb94ee6SHong Zhang /* construct I-gamma*Jac */ 607cb94ee6SHong Zhang gm = -_gamma; 617cb94ee6SHong Zhang ierr = MatScale(Jac,gm);CHKERRQ(ierr); 627cb94ee6SHong Zhang ierr = MatShift(Jac,one);CHKERRQ(ierr); 637cb94ee6SHong Zhang 647cb94ee6SHong Zhang ierr = PCSetOperators(pc,Jac,Jac,str);CHKERRQ(ierr); 657cb94ee6SHong Zhang PetscFunctionReturn(0); 667cb94ee6SHong Zhang } 677cb94ee6SHong Zhang 687cb94ee6SHong Zhang /* 697cb94ee6SHong Zhang TSPSolve_Sundials - routine that we provide to Sundials that applies the preconditioner. 707cb94ee6SHong Zhang */ 717cb94ee6SHong Zhang #undef __FUNCT__ 727cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials" 734ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z, 744ac7836bSHong Zhang realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp) 757cb94ee6SHong Zhang { 767cb94ee6SHong Zhang TS ts = (TS) P_data; 777cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 787cb94ee6SHong Zhang PC pc = cvode->pc; 797cb94ee6SHong Zhang Vec rr = cvode->w1,zz = cvode->w2; 807cb94ee6SHong Zhang PetscErrorCode ierr; 817cb94ee6SHong Zhang PetscScalar *r_data,*z_data; 827cb94ee6SHong Zhang 837cb94ee6SHong Zhang PetscFunctionBegin; 847cb94ee6SHong Zhang /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/ 857cb94ee6SHong Zhang r_data = (PetscScalar *) N_VGetArrayPointer(r); 867cb94ee6SHong Zhang z_data = (PetscScalar *) N_VGetArrayPointer(z); 877cb94ee6SHong Zhang ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr); 887cb94ee6SHong Zhang ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr); 894ac7836bSHong Zhang 907cb94ee6SHong Zhang /* Solve the Px=r and put the result in zz */ 917cb94ee6SHong Zhang ierr = PCApply(pc,rr,zz); CHKERRQ(ierr); 927cb94ee6SHong Zhang ierr = VecResetArray(rr); CHKERRQ(ierr); 937cb94ee6SHong Zhang ierr = VecResetArray(zz); CHKERRQ(ierr); 947cb94ee6SHong Zhang PetscFunctionReturn(0); 957cb94ee6SHong Zhang } 967cb94ee6SHong Zhang 977cb94ee6SHong Zhang /* 987cb94ee6SHong Zhang TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side. 997cb94ee6SHong Zhang */ 1007cb94ee6SHong Zhang #undef __FUNCT__ 1017cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials" 1024ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx) 1037cb94ee6SHong Zhang { 1047cb94ee6SHong Zhang TS ts = (TS) ctx; 1057adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)ts)->comm; 1067cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1077cb94ee6SHong Zhang Vec yy = cvode->w1,yyd = cvode->w2; 1087cb94ee6SHong Zhang PetscScalar *y_data,*ydot_data; 1097cb94ee6SHong Zhang PetscErrorCode ierr; 1107cb94ee6SHong Zhang 1117cb94ee6SHong Zhang PetscFunctionBegin; 1125ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 1137cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(y); 1147cb94ee6SHong Zhang ydot_data = (PetscScalar *) N_VGetArrayPointer(ydot); 1154ac7836bSHong Zhang ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr); 1164ac7836bSHong Zhang ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr); 1174ac7836bSHong Zhang 1187cb94ee6SHong Zhang /* now compute the right hand side function */ 1197cb94ee6SHong Zhang ierr = TSComputeRHSFunction(ts,t,yy,yyd); CHKERRABORT(comm,ierr); 1207cb94ee6SHong Zhang ierr = VecResetArray(yy); CHKERRABORT(comm,ierr); 1217cb94ee6SHong Zhang ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr); 1224ac7836bSHong Zhang PetscFunctionReturn(0); 1237cb94ee6SHong Zhang } 1247cb94ee6SHong Zhang 1257cb94ee6SHong Zhang /* 1267cb94ee6SHong Zhang TSStep_Sundials_Nonlinear - Calls Sundials to integrate the ODE. 1277cb94ee6SHong Zhang */ 1287cb94ee6SHong Zhang #undef __FUNCT__ 1297cb94ee6SHong Zhang #define __FUNCT__ "TSStep_Sundials_Nonlinear" 1307cb94ee6SHong Zhang PetscErrorCode TSStep_Sundials_Nonlinear(TS ts,int *steps,double *time) 1317cb94ee6SHong Zhang { 1327cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1337cb94ee6SHong Zhang Vec sol = ts->vec_sol; 1347cb94ee6SHong Zhang PetscErrorCode ierr; 1355d47aee6SHong Zhang PetscInt i,max_steps = ts->max_steps,flag; 1369f94935aSHong Zhang long int its,nsteps; 1377cb94ee6SHong Zhang realtype t,tout; 1387cb94ee6SHong Zhang PetscScalar *y_data; 1397cb94ee6SHong Zhang void *mem; 1407cb94ee6SHong Zhang 1417cb94ee6SHong Zhang PetscFunctionBegin; 14216016685SHong Zhang mem = cvode->mem; 1437cb94ee6SHong Zhang tout = ts->max_time; 1447cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr); 1457cb94ee6SHong Zhang N_VSetArrayPointer((realtype *)y_data,cvode->y); 1467cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 1477cb94ee6SHong Zhang for (i = 0; i < max_steps; i++) { 1487cb94ee6SHong Zhang if (ts->ptime >= ts->max_time) break; 1493f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 1502bfc04deSHong Zhang if (cvode->monitorstep){ 15128adb3f7SHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP); 1522bfc04deSHong Zhang } else { 1532bfc04deSHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL); 1542bfc04deSHong Zhang } 1559f94935aSHong Zhang 1569f94935aSHong Zhang if (flag){ /* display error message */ 1579f94935aSHong Zhang switch (flag){ 1589f94935aSHong Zhang case CV_ILL_INPUT: 1599f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT"); 1609f94935aSHong Zhang break; 1619f94935aSHong Zhang case CV_TOO_CLOSE: 1629f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE"); 1639f94935aSHong Zhang break; 1643c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1659f94935aSHong Zhang PetscReal tcur; 1669f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 1679f94935aSHong Zhang ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr); 1689f94935aSHong Zhang SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_WORK. At t=%G, nsteps %D exceeds mxstep %D. Increase '-ts_max_steps <>' or modify TSSetDuration()",tcur,nsteps,ts->max_steps); 1693c7fefeeSJed Brown } break; 1709f94935aSHong Zhang case CV_TOO_MUCH_ACC: 1719f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC"); 1729f94935aSHong Zhang break; 1739f94935aSHong Zhang case CV_ERR_FAILURE: 1749f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE"); 1759f94935aSHong Zhang break; 1769f94935aSHong Zhang case CV_CONV_FAILURE: 1779f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE"); 1789f94935aSHong Zhang break; 1799f94935aSHong Zhang case CV_LINIT_FAIL: 1809f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL"); 1819f94935aSHong Zhang break; 1829f94935aSHong Zhang case CV_LSETUP_FAIL: 1839f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL"); 1849f94935aSHong Zhang break; 1859f94935aSHong Zhang case CV_LSOLVE_FAIL: 1869f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL"); 1879f94935aSHong Zhang break; 1889f94935aSHong Zhang case CV_RHSFUNC_FAIL: 1899f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL"); 1909f94935aSHong Zhang break; 1919f94935aSHong Zhang case CV_FIRST_RHSFUNC_ERR: 1929f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR"); 1939f94935aSHong Zhang break; 1949f94935aSHong Zhang case CV_REPTD_RHSFUNC_ERR: 1959f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR"); 1969f94935aSHong Zhang break; 1979f94935aSHong Zhang case CV_UNREC_RHSFUNC_ERR: 1989f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR"); 1999f94935aSHong Zhang break; 2009f94935aSHong Zhang case CV_RTFUNC_FAIL: 2019f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL"); 2029f94935aSHong Zhang break; 2039f94935aSHong Zhang default: 2049f94935aSHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag); 2059f94935aSHong Zhang } 2069f94935aSHong Zhang } 2079f94935aSHong Zhang 2087cb94ee6SHong Zhang if (t > ts->max_time && cvode->exact_final_time) { 2097cb94ee6SHong Zhang /* interpolate to final requested time */ 2107cb94ee6SHong Zhang ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr); 2117cb94ee6SHong Zhang t = tout; 2127cb94ee6SHong Zhang } 2137cb94ee6SHong Zhang ts->time_step = t - ts->ptime; 2147cb94ee6SHong Zhang ts->ptime = t; 2157cb94ee6SHong Zhang 2167cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 2177cb94ee6SHong Zhang ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr); 2187cb94ee6SHong Zhang ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr); 2197cb94ee6SHong Zhang ierr = VecResetArray(cvode->w1); CHKERRQ(ierr); 2207cb94ee6SHong Zhang ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr); 2217cb94ee6SHong Zhang ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr); 2227cb94ee6SHong Zhang ts->nonlinear_its = its; 2234ac7836bSHong Zhang ierr = CVSpilsGetNumLinIters(mem, &its); 2247cb94ee6SHong Zhang ts->linear_its = its; 2257cb94ee6SHong Zhang ts->steps++; 2263f2090d5SJed Brown ierr = TSPostStep(ts);CHKERRQ(ierr); 2277cb94ee6SHong Zhang ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr); 2287cb94ee6SHong Zhang } 2299f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 2309f94935aSHong Zhang *steps = nsteps; 2317cb94ee6SHong Zhang *time = t; 2327cb94ee6SHong Zhang PetscFunctionReturn(0); 2337cb94ee6SHong Zhang } 2347cb94ee6SHong Zhang 2357cb94ee6SHong Zhang #undef __FUNCT__ 236*277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Sundials" 237*277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts) 238*277b19d0SLisandro Dalcin { 239*277b19d0SLisandro Dalcin TS_Sundials *cvode = (TS_Sundials*)ts->data; 240*277b19d0SLisandro Dalcin PetscErrorCode ierr; 241*277b19d0SLisandro Dalcin 242*277b19d0SLisandro Dalcin PetscFunctionBegin; 243*277b19d0SLisandro Dalcin if (cvode->pc) {ierr = PCReset(cvode->pc);CHKERRQ(ierr);} 244*277b19d0SLisandro Dalcin if (cvode->pmat) {ierr = MatDestroy(cvode->pmat);CHKERRQ(ierr);} 245*277b19d0SLisandro Dalcin if (cvode->update) {ierr = VecDestroy(cvode->update);CHKERRQ(ierr);} 246*277b19d0SLisandro Dalcin if (cvode->func) {ierr = VecDestroy(cvode->func);CHKERRQ(ierr);} 247*277b19d0SLisandro Dalcin if (cvode->rhs) {ierr = VecDestroy(cvode->rhs);CHKERRQ(ierr);} 248*277b19d0SLisandro Dalcin if (cvode->w1) {ierr = VecDestroy(cvode->w1);CHKERRQ(ierr);} 249*277b19d0SLisandro Dalcin if (cvode->w2) {ierr = VecDestroy(cvode->w2);CHKERRQ(ierr);} 250*277b19d0SLisandro Dalcin if (cvode->mem) {CVodeFree(&cvode->mem);} 251*277b19d0SLisandro Dalcin PetscFunctionReturn(0); 252*277b19d0SLisandro Dalcin } 253*277b19d0SLisandro Dalcin 254*277b19d0SLisandro Dalcin #undef __FUNCT__ 2557cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials" 2567cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts) 2577cb94ee6SHong Zhang { 2587cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2597cb94ee6SHong Zhang PetscErrorCode ierr; 2607cb94ee6SHong Zhang 2617cb94ee6SHong Zhang PetscFunctionBegin; 262*277b19d0SLisandro Dalcin ierr = TSReset_Sundials(ts);CHKERRQ(ierr); 2637cb94ee6SHong Zhang if (cvode->pc) {ierr = PCDestroy(cvode->pc);CHKERRQ(ierr);} 264a07356b8SSatish Balay ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr); 265*277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 2667cb94ee6SHong Zhang PetscFunctionReturn(0); 2677cb94ee6SHong Zhang } 2687cb94ee6SHong Zhang 2697cb94ee6SHong Zhang #undef __FUNCT__ 2707cb94ee6SHong Zhang #define __FUNCT__ "TSSetUp_Sundials_Nonlinear" 2717cb94ee6SHong Zhang PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts) 2727cb94ee6SHong Zhang { 2737cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2747cb94ee6SHong Zhang PetscErrorCode ierr; 27516016685SHong Zhang PetscInt glosize,locsize,i,flag; 2767cb94ee6SHong Zhang PetscScalar *y_data,*parray; 27716016685SHong Zhang void *mem; 27816016685SHong Zhang const PCType pctype; 279ace3abfcSBarry Smith PetscBool pcnone; 28016016685SHong Zhang Vec sol = ts->vec_sol; 2817cb94ee6SHong Zhang 2827cb94ee6SHong Zhang PetscFunctionBegin; 2837cb94ee6SHong Zhang ierr = PCSetFromOptions(cvode->pc);CHKERRQ(ierr); 2849f94935aSHong Zhang 2857cb94ee6SHong Zhang /* get the vector size */ 2867cb94ee6SHong Zhang ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr); 2877cb94ee6SHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr); 2887cb94ee6SHong Zhang 2897cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 290a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 291e32f2f54SBarry Smith if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated"); 2927cb94ee6SHong Zhang 29328adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 2947cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr); 2957cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y); 2967cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 2977cb94ee6SHong Zhang /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/ 2987cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 2997cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr); 3007cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->func);CHKERRQ(ierr); 3017cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr); 3027cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->func);CHKERRQ(ierr); 3037cb94ee6SHong Zhang 3047cb94ee6SHong Zhang /* 3057cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 3067cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 3077cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 3087cb94ee6SHong Zhang */ 3097adad957SLisandro Dalcin ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr); 3107adad957SLisandro Dalcin ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr); 3117cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr); 3127cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr); 31316016685SHong Zhang 31416016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 31516016685SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 316e32f2f54SBarry Smith if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails"); 31716016685SHong Zhang cvode->mem = mem; 31816016685SHong Zhang 31916016685SHong Zhang /* Set the pointer to user-defined data */ 32016016685SHong Zhang flag = CVodeSetUserData(mem, ts); 321e32f2f54SBarry Smith if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails"); 32216016685SHong Zhang 323fc6b9e64SJed Brown /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 324fc6b9e64SJed Brown flag = CVodeSetInitStep(mem,(realtype)ts->initial_time_step); 325fc6b9e64SJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed"); 326f1cd61daSJed Brown if (cvode->mindt > 0) { 327f1cd61daSJed Brown flag = CVodeSetMinStep(mem,(realtype)cvode->mindt); 3289f94935aSHong Zhang if (flag){ 3299f94935aSHong Zhang if (flag == CV_MEM_NULL){ 3309f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 3319f94935aSHong Zhang } else if (flag == CV_ILL_INPUT){ 3329f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size"); 3339f94935aSHong Zhang } else { 3349f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed"); 3359f94935aSHong Zhang } 3369f94935aSHong Zhang } 337f1cd61daSJed Brown } 338f1cd61daSJed Brown if (cvode->maxdt > 0) { 339f1cd61daSJed Brown flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt); 340f1cd61daSJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed"); 341f1cd61daSJed Brown } 342f1cd61daSJed Brown 34316016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 34416016685SHong Zhang * user's right hand side function in u'=f(t,u), the inital time T0, and 34516016685SHong Zhang * the initial dependent variable vector cvode->y */ 34616016685SHong Zhang flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 34716016685SHong Zhang if (flag){ 348e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag); 34916016685SHong Zhang } 35016016685SHong Zhang 3519f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 35216016685SHong Zhang flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 35316016685SHong Zhang if (flag){ 354e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag); 35516016685SHong Zhang } 35616016685SHong Zhang 3579f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 3589f94935aSHong Zhang flag = CVodeSetMaxNumSteps(mem,ts->max_steps); 35916016685SHong Zhang ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 36016016685SHong Zhang 36116016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 36216016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 36316016685SHong Zhang ierr = PCGetType(cvode->pc,&pctype);CHKERRQ(ierr); 36416016685SHong Zhang ierr = PetscTypeCompare((PetscObject)cvode->pc,PCNONE,&pcnone);CHKERRQ(ierr); 36516016685SHong Zhang if (pcnone){ 36616016685SHong Zhang flag = CVSpgmr(mem,PREC_NONE,0); 367e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 36816016685SHong Zhang } else { 36916016685SHong Zhang flag = CVSpgmr(mem,PREC_LEFT,0); 370e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 37116016685SHong Zhang 37216016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 37316016685SHong Zhang and the pointer to the user-defined block data */ 37416016685SHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 375e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 37616016685SHong Zhang } 37716016685SHong Zhang 37816016685SHong Zhang flag = CVSpilsSetGSType(mem, MODIFIED_GS); 379e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag); 3807cb94ee6SHong Zhang PetscFunctionReturn(0); 3817cb94ee6SHong Zhang } 3827cb94ee6SHong Zhang 3836fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 384dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0}; 3856fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 386dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0}; 387a04cf4d8SBarry Smith 3887cb94ee6SHong Zhang #undef __FUNCT__ 3897cb94ee6SHong Zhang #define __FUNCT__ "TSSetFromOptions_Sundials_Nonlinear" 3907cb94ee6SHong Zhang PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts) 3917cb94ee6SHong Zhang { 3927cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 3937cb94ee6SHong Zhang PetscErrorCode ierr; 3947cb94ee6SHong Zhang int indx; 395ace3abfcSBarry Smith PetscBool flag; 3967cb94ee6SHong Zhang 3977cb94ee6SHong Zhang PetscFunctionBegin; 3987cb94ee6SHong Zhang ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr); 3996fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr); 4007cb94ee6SHong Zhang if (flag) { 4016fadb2cdSHong Zhang ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr); 4027cb94ee6SHong Zhang } 4036fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr); 4047cb94ee6SHong Zhang if (flag) { 4057cb94ee6SHong Zhang ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 4067cb94ee6SHong Zhang } 4077cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr); 4087cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr); 409f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr); 410f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr); 4117cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr); 4127cb94ee6SHong Zhang ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr); 413acfcf0e5SJed Brown ierr = PetscOptionsBool("-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); 414acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr); 4157cb94ee6SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4167cb94ee6SHong Zhang PetscFunctionReturn(0); 4177cb94ee6SHong Zhang } 4187cb94ee6SHong Zhang 4197cb94ee6SHong Zhang #undef __FUNCT__ 4207cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials" 4217cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 4227cb94ee6SHong Zhang { 4237cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4247cb94ee6SHong Zhang PetscErrorCode ierr; 4257cb94ee6SHong Zhang char *type; 4267cb94ee6SHong Zhang char atype[] = "Adams"; 4277cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 428ace3abfcSBarry Smith PetscBool iascii,isstring; 4292c823083SHong Zhang long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 4302c823083SHong Zhang PetscInt qlast,qcur; 4315d47aee6SHong Zhang PetscReal hinused,hlast,hcur,tcur,tolsfac; 4327cb94ee6SHong Zhang 4337cb94ee6SHong Zhang PetscFunctionBegin; 4347cb94ee6SHong Zhang if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;} 4357cb94ee6SHong Zhang else {type = btype;} 4367cb94ee6SHong Zhang 4372692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4382692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 4397cb94ee6SHong Zhang if (iascii) { 4407cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr); 4417cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr); 4427cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr); 4437cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr); 4447cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr); 4457cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 4467cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4477cb94ee6SHong Zhang } else { 4487cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4497cb94ee6SHong Zhang } 450f1cd61daSJed Brown if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);} 451f1cd61daSJed Brown if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);} 4522c823083SHong Zhang 4535d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4545d47aee6SHong Zhang ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr); 4555d47aee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr); 4562c823083SHong Zhang ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals, 4572c823083SHong Zhang &nlinsetups,&nfails,&qlast,&qcur, 4582c823083SHong Zhang &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr); 4592c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr); 4602c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr); 4612c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr); 4622c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr); 4632c823083SHong Zhang 4642c823083SHong Zhang ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr); 4652c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr); 4662c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr); 4672c823083SHong Zhang 4682c823083SHong Zhang ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */ 4692c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr); 4702c823083SHong Zhang ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr); 4712c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr); 4722c823083SHong Zhang ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4732c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr); 4742c823083SHong Zhang ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr); 4752c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr); 4762c823083SHong Zhang ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4772c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr); 4782c823083SHong Zhang ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4792c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr); 4807cb94ee6SHong Zhang } else if (isstring) { 4817cb94ee6SHong Zhang ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 4827cb94ee6SHong Zhang } else { 483e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name); 4847cb94ee6SHong Zhang } 4857cb94ee6SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 4867cb94ee6SHong Zhang ierr = PCView(cvode->pc,viewer);CHKERRQ(ierr); 4877cb94ee6SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 4887cb94ee6SHong Zhang PetscFunctionReturn(0); 4897cb94ee6SHong Zhang } 4907cb94ee6SHong Zhang 4917cb94ee6SHong Zhang 4927cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 4937cb94ee6SHong Zhang EXTERN_C_BEGIN 4947cb94ee6SHong Zhang #undef __FUNCT__ 4957cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials" 4967087cfbeSBarry Smith PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 4977cb94ee6SHong Zhang { 4987cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4997cb94ee6SHong Zhang 5007cb94ee6SHong Zhang PetscFunctionBegin; 5017cb94ee6SHong Zhang cvode->cvode_type = type; 5027cb94ee6SHong Zhang PetscFunctionReturn(0); 5037cb94ee6SHong Zhang } 5047cb94ee6SHong Zhang EXTERN_C_END 5057cb94ee6SHong Zhang 5067cb94ee6SHong Zhang EXTERN_C_BEGIN 5077cb94ee6SHong Zhang #undef __FUNCT__ 5087cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials" 5097087cfbeSBarry Smith PetscErrorCode TSSundialsSetGMRESRestart_Sundials(TS ts,int restart) 5107cb94ee6SHong Zhang { 5117cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5127cb94ee6SHong Zhang 5137cb94ee6SHong Zhang PetscFunctionBegin; 5147cb94ee6SHong Zhang cvode->restart = restart; 5157cb94ee6SHong Zhang PetscFunctionReturn(0); 5167cb94ee6SHong Zhang } 5177cb94ee6SHong Zhang EXTERN_C_END 5187cb94ee6SHong Zhang 5197cb94ee6SHong Zhang EXTERN_C_BEGIN 5207cb94ee6SHong Zhang #undef __FUNCT__ 5217cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials" 5227087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,double tol) 5237cb94ee6SHong Zhang { 5247cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5257cb94ee6SHong Zhang 5267cb94ee6SHong Zhang PetscFunctionBegin; 5277cb94ee6SHong Zhang cvode->linear_tol = tol; 5287cb94ee6SHong Zhang PetscFunctionReturn(0); 5297cb94ee6SHong Zhang } 5307cb94ee6SHong Zhang EXTERN_C_END 5317cb94ee6SHong Zhang 5327cb94ee6SHong Zhang EXTERN_C_BEGIN 5337cb94ee6SHong Zhang #undef __FUNCT__ 5347cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials" 5357087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 5367cb94ee6SHong Zhang { 5377cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5387cb94ee6SHong Zhang 5397cb94ee6SHong Zhang PetscFunctionBegin; 5407cb94ee6SHong Zhang cvode->gtype = type; 5417cb94ee6SHong Zhang PetscFunctionReturn(0); 5427cb94ee6SHong Zhang } 5437cb94ee6SHong Zhang EXTERN_C_END 5447cb94ee6SHong Zhang 5457cb94ee6SHong Zhang EXTERN_C_BEGIN 5467cb94ee6SHong Zhang #undef __FUNCT__ 5477cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials" 5487087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel) 5497cb94ee6SHong Zhang { 5507cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5517cb94ee6SHong Zhang 5527cb94ee6SHong Zhang PetscFunctionBegin; 5537cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 5547cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 5557cb94ee6SHong Zhang PetscFunctionReturn(0); 5567cb94ee6SHong Zhang } 5577cb94ee6SHong Zhang EXTERN_C_END 5587cb94ee6SHong Zhang 5597cb94ee6SHong Zhang EXTERN_C_BEGIN 5607cb94ee6SHong Zhang #undef __FUNCT__ 561f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials" 5627087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 563f1cd61daSJed Brown { 564f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 565f1cd61daSJed Brown 566f1cd61daSJed Brown PetscFunctionBegin; 567f1cd61daSJed Brown cvode->mindt = mindt; 568f1cd61daSJed Brown PetscFunctionReturn(0); 569f1cd61daSJed Brown } 570f1cd61daSJed Brown EXTERN_C_END 571f1cd61daSJed Brown 572f1cd61daSJed Brown EXTERN_C_BEGIN 573f1cd61daSJed Brown #undef __FUNCT__ 574f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials" 5757087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 576f1cd61daSJed Brown { 577f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 578f1cd61daSJed Brown 579f1cd61daSJed Brown PetscFunctionBegin; 580f1cd61daSJed Brown cvode->maxdt = maxdt; 581f1cd61daSJed Brown PetscFunctionReturn(0); 582f1cd61daSJed Brown } 583f1cd61daSJed Brown EXTERN_C_END 584f1cd61daSJed Brown 585f1cd61daSJed Brown EXTERN_C_BEGIN 586f1cd61daSJed Brown #undef __FUNCT__ 5877cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials" 5887087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc) 5897cb94ee6SHong Zhang { 5907cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5917cb94ee6SHong Zhang 5927cb94ee6SHong Zhang PetscFunctionBegin; 5937cb94ee6SHong Zhang *pc = cvode->pc; 5947cb94ee6SHong Zhang PetscFunctionReturn(0); 5957cb94ee6SHong Zhang } 5967cb94ee6SHong Zhang EXTERN_C_END 5977cb94ee6SHong Zhang 5987cb94ee6SHong Zhang EXTERN_C_BEGIN 5997cb94ee6SHong Zhang #undef __FUNCT__ 6007cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials" 6017087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 6027cb94ee6SHong Zhang { 6037cb94ee6SHong Zhang PetscFunctionBegin; 6042c823083SHong Zhang if (nonlin) *nonlin = ts->nonlinear_its; 6052c823083SHong Zhang if (lin) *lin = ts->linear_its; 6067cb94ee6SHong Zhang PetscFunctionReturn(0); 6077cb94ee6SHong Zhang } 6087cb94ee6SHong Zhang EXTERN_C_END 6097cb94ee6SHong Zhang 6107cb94ee6SHong Zhang EXTERN_C_BEGIN 6117cb94ee6SHong Zhang #undef __FUNCT__ 6127cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials" 6137087cfbeSBarry Smith PetscErrorCode TSSundialsSetExactFinalTime_Sundials(TS ts,PetscBool s) 6147cb94ee6SHong Zhang { 6157cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 6167cb94ee6SHong Zhang 6177cb94ee6SHong Zhang PetscFunctionBegin; 6187cb94ee6SHong Zhang cvode->exact_final_time = s; 6197cb94ee6SHong Zhang PetscFunctionReturn(0); 6207cb94ee6SHong Zhang } 6217cb94ee6SHong Zhang EXTERN_C_END 6222bfc04deSHong Zhang 6232bfc04deSHong Zhang EXTERN_C_BEGIN 6242bfc04deSHong Zhang #undef __FUNCT__ 6252bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials" 6267087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s) 6272bfc04deSHong Zhang { 6282bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 6292bfc04deSHong Zhang 6302bfc04deSHong Zhang PetscFunctionBegin; 6312bfc04deSHong Zhang cvode->monitorstep = s; 6322bfc04deSHong Zhang PetscFunctionReturn(0); 6332bfc04deSHong Zhang } 6342bfc04deSHong Zhang EXTERN_C_END 6357cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 6367cb94ee6SHong Zhang 6377cb94ee6SHong Zhang #undef __FUNCT__ 6387cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations" 6397cb94ee6SHong Zhang /*@C 6407cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 6417cb94ee6SHong Zhang 6427cb94ee6SHong Zhang Not Collective 6437cb94ee6SHong Zhang 6447cb94ee6SHong Zhang Input parameters: 6457cb94ee6SHong Zhang . ts - the time-step context 6467cb94ee6SHong Zhang 6477cb94ee6SHong Zhang Output Parameters: 6487cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 6497cb94ee6SHong Zhang - lin - number of linear iterations 6507cb94ee6SHong Zhang 6517cb94ee6SHong Zhang Level: advanced 6527cb94ee6SHong Zhang 6537cb94ee6SHong Zhang Notes: 6547cb94ee6SHong Zhang These return the number since the creation of the TS object 6557cb94ee6SHong Zhang 6567cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations 6577cb94ee6SHong Zhang 6587cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6597cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 6607cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 661f1cd61daSJed Brown TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSundialsSetExactFinalTime() 6627cb94ee6SHong Zhang 6637cb94ee6SHong Zhang @*/ 6647087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 6657cb94ee6SHong Zhang { 6664ac538c5SBarry Smith PetscErrorCode ierr; 6677cb94ee6SHong Zhang 6687cb94ee6SHong Zhang PetscFunctionBegin; 6694ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr); 6707cb94ee6SHong Zhang PetscFunctionReturn(0); 6717cb94ee6SHong Zhang } 6727cb94ee6SHong Zhang 6737cb94ee6SHong Zhang #undef __FUNCT__ 6747cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType" 6757cb94ee6SHong Zhang /*@ 6767cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 6777cb94ee6SHong Zhang 6783f9fe445SBarry Smith Logically Collective on TS 6797cb94ee6SHong Zhang 6807cb94ee6SHong Zhang Input parameters: 6817cb94ee6SHong Zhang + ts - the time-step context 6827cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 6837cb94ee6SHong Zhang 6847cb94ee6SHong Zhang Level: intermediate 6857cb94ee6SHong Zhang 6867cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula 6877cb94ee6SHong Zhang 6887cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetGMRESRestart(), 6897cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 6907cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6917cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 6927cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 6937cb94ee6SHong Zhang @*/ 6947087cfbeSBarry Smith PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type) 6957cb94ee6SHong Zhang { 6964ac538c5SBarry Smith PetscErrorCode ierr; 6977cb94ee6SHong Zhang 6987cb94ee6SHong Zhang PetscFunctionBegin; 6994ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr); 7007cb94ee6SHong Zhang PetscFunctionReturn(0); 7017cb94ee6SHong Zhang } 7027cb94ee6SHong Zhang 7037cb94ee6SHong Zhang #undef __FUNCT__ 7047cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart" 7057cb94ee6SHong Zhang /*@ 7067cb94ee6SHong Zhang TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by 7077cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 7087cb94ee6SHong Zhang this is ALSO the maximum number of GMRES steps that will be used. 7097cb94ee6SHong Zhang 7103f9fe445SBarry Smith Logically Collective on TS 7117cb94ee6SHong Zhang 7127cb94ee6SHong Zhang Input parameters: 7137cb94ee6SHong Zhang + ts - the time-step context 7147cb94ee6SHong Zhang - restart - number of direction vectors (the restart size). 7157cb94ee6SHong Zhang 7167cb94ee6SHong Zhang Level: advanced 7177cb94ee6SHong Zhang 7187cb94ee6SHong Zhang .keywords: GMRES, restart 7197cb94ee6SHong Zhang 7207cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 7217cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 7227cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7237cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 7247cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 7257cb94ee6SHong Zhang 7267cb94ee6SHong Zhang @*/ 7277087cfbeSBarry Smith PetscErrorCode TSSundialsSetGMRESRestart(TS ts,int restart) 7287cb94ee6SHong Zhang { 7294ac538c5SBarry Smith PetscErrorCode ierr; 7307cb94ee6SHong Zhang 7317cb94ee6SHong Zhang PetscFunctionBegin; 732c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(ts,restart,2); 7334ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetGMRESRestart_C",(TS,int),(ts,restart));CHKERRQ(ierr); 7347cb94ee6SHong Zhang PetscFunctionReturn(0); 7357cb94ee6SHong Zhang } 7367cb94ee6SHong Zhang 7377cb94ee6SHong Zhang #undef __FUNCT__ 7387cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance" 7397cb94ee6SHong Zhang /*@ 7407cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 7417cb94ee6SHong Zhang system by SUNDIALS. 7427cb94ee6SHong Zhang 7433f9fe445SBarry Smith Logically Collective on TS 7447cb94ee6SHong Zhang 7457cb94ee6SHong Zhang Input parameters: 7467cb94ee6SHong Zhang + ts - the time-step context 7477cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 7487cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 7497cb94ee6SHong Zhang 7507cb94ee6SHong Zhang Level: advanced 7517cb94ee6SHong Zhang 7527cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS 7537cb94ee6SHong Zhang 7547cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7557cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 7567cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7577cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 7587cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 7597cb94ee6SHong Zhang 7607cb94ee6SHong Zhang @*/ 7617087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance(TS ts,double tol) 7627cb94ee6SHong Zhang { 7634ac538c5SBarry Smith PetscErrorCode ierr; 7647cb94ee6SHong Zhang 7657cb94ee6SHong Zhang PetscFunctionBegin; 766c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,tol,2); 7674ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr); 7687cb94ee6SHong Zhang PetscFunctionReturn(0); 7697cb94ee6SHong Zhang } 7707cb94ee6SHong Zhang 7717cb94ee6SHong Zhang #undef __FUNCT__ 7727cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType" 7737cb94ee6SHong Zhang /*@ 7747cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 7757cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 7767cb94ee6SHong Zhang 7773f9fe445SBarry Smith Logically Collective on TS 7787cb94ee6SHong Zhang 7797cb94ee6SHong Zhang Input parameters: 7807cb94ee6SHong Zhang + ts - the time-step context 7817cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 7827cb94ee6SHong Zhang 7837cb94ee6SHong Zhang Level: advanced 7847cb94ee6SHong Zhang 7857cb94ee6SHong Zhang .keywords: Sundials, orthogonalization 7867cb94ee6SHong Zhang 7877cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7887cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 7897cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7907cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 7917cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 7927cb94ee6SHong Zhang 7937cb94ee6SHong Zhang @*/ 7947087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 7957cb94ee6SHong Zhang { 7964ac538c5SBarry Smith PetscErrorCode ierr; 7977cb94ee6SHong Zhang 7987cb94ee6SHong Zhang PetscFunctionBegin; 7994ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr); 8007cb94ee6SHong Zhang PetscFunctionReturn(0); 8017cb94ee6SHong Zhang } 8027cb94ee6SHong Zhang 8037cb94ee6SHong Zhang #undef __FUNCT__ 8047cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance" 8057cb94ee6SHong Zhang /*@ 8067cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 8077cb94ee6SHong Zhang Sundials for error control. 8087cb94ee6SHong Zhang 8093f9fe445SBarry Smith Logically Collective on TS 8107cb94ee6SHong Zhang 8117cb94ee6SHong Zhang Input parameters: 8127cb94ee6SHong Zhang + ts - the time-step context 8137cb94ee6SHong Zhang . aabs - the absolute tolerance 8147cb94ee6SHong Zhang - rel - the relative tolerance 8157cb94ee6SHong Zhang 8167cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 8177cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 8187cb94ee6SHong Zhang 8197cb94ee6SHong Zhang Level: intermediate 8207cb94ee6SHong Zhang 8217cb94ee6SHong Zhang .keywords: Sundials, tolerance 8227cb94ee6SHong Zhang 8237cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8247cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 8257cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8267cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 8277cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 8287cb94ee6SHong Zhang 8297cb94ee6SHong Zhang @*/ 8307087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance(TS ts,double aabs,double rel) 8317cb94ee6SHong Zhang { 8324ac538c5SBarry Smith PetscErrorCode ierr; 8337cb94ee6SHong Zhang 8347cb94ee6SHong Zhang PetscFunctionBegin; 8354ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr); 8367cb94ee6SHong Zhang PetscFunctionReturn(0); 8377cb94ee6SHong Zhang } 8387cb94ee6SHong Zhang 8397cb94ee6SHong Zhang #undef __FUNCT__ 8407cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC" 8417cb94ee6SHong Zhang /*@ 8427cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 8437cb94ee6SHong Zhang 8447cb94ee6SHong Zhang Input Parameter: 8457cb94ee6SHong Zhang . ts - the time-step context 8467cb94ee6SHong Zhang 8477cb94ee6SHong Zhang Output Parameter: 8487cb94ee6SHong Zhang . pc - the preconditioner context 8497cb94ee6SHong Zhang 8507cb94ee6SHong Zhang Level: advanced 8517cb94ee6SHong Zhang 8527cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8537cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 8547cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8557cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 8567cb94ee6SHong Zhang @*/ 8577087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC(TS ts,PC *pc) 8587cb94ee6SHong Zhang { 8594ac538c5SBarry Smith PetscErrorCode ierr; 8607cb94ee6SHong Zhang 8617cb94ee6SHong Zhang PetscFunctionBegin; 8624ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr); 8637cb94ee6SHong Zhang PetscFunctionReturn(0); 8647cb94ee6SHong Zhang } 8657cb94ee6SHong Zhang 8667cb94ee6SHong Zhang #undef __FUNCT__ 8677cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime" 8687cb94ee6SHong Zhang /*@ 8697cb94ee6SHong Zhang TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the 8707cb94ee6SHong Zhang exact final time requested by the user or just returns it at the final time 8717cb94ee6SHong Zhang it computed. (Defaults to true). 8727cb94ee6SHong Zhang 8737cb94ee6SHong Zhang Input Parameter: 8747cb94ee6SHong Zhang + ts - the time-step context 8757cb94ee6SHong Zhang - ft - PETSC_TRUE if interpolates, else PETSC_FALSE 8767cb94ee6SHong Zhang 8777cb94ee6SHong Zhang Level: beginner 8787cb94ee6SHong Zhang 8797cb94ee6SHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8807cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 8817cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8827cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 8837cb94ee6SHong Zhang @*/ 8847087cfbeSBarry Smith PetscErrorCode TSSundialsSetExactFinalTime(TS ts,PetscBool ft) 8857cb94ee6SHong Zhang { 8864ac538c5SBarry Smith PetscErrorCode ierr; 8877cb94ee6SHong Zhang 8887cb94ee6SHong Zhang PetscFunctionBegin; 8894ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetExactFinalTime_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 8907cb94ee6SHong Zhang PetscFunctionReturn(0); 8917cb94ee6SHong Zhang } 8927cb94ee6SHong Zhang 8932bfc04deSHong Zhang #undef __FUNCT__ 894f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep" 895f1cd61daSJed Brown /*@ 896f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 897f1cd61daSJed Brown 898f1cd61daSJed Brown Input Parameter: 899f1cd61daSJed Brown + ts - the time-step context 900f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 901f1cd61daSJed Brown 902fc6b9e64SJed Brown Note: 903fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 904fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 905fc6b9e64SJed Brown 906f1cd61daSJed Brown Level: beginner 907f1cd61daSJed Brown 908f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 909f1cd61daSJed Brown @*/ 9107087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 911f1cd61daSJed Brown { 9124ac538c5SBarry Smith PetscErrorCode ierr; 913f1cd61daSJed Brown 914f1cd61daSJed Brown PetscFunctionBegin; 9154ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr); 916f1cd61daSJed Brown PetscFunctionReturn(0); 917f1cd61daSJed Brown } 918f1cd61daSJed Brown 919f1cd61daSJed Brown #undef __FUNCT__ 920f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep" 921f1cd61daSJed Brown /*@ 922f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 923f1cd61daSJed Brown 924f1cd61daSJed Brown Input Parameter: 925f1cd61daSJed Brown + ts - the time-step context 926f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 927f1cd61daSJed Brown 928f1cd61daSJed Brown Level: beginner 929f1cd61daSJed Brown 930f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 931f1cd61daSJed Brown @*/ 9327087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 933f1cd61daSJed Brown { 9344ac538c5SBarry Smith PetscErrorCode ierr; 935f1cd61daSJed Brown 936f1cd61daSJed Brown PetscFunctionBegin; 9374ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 938f1cd61daSJed Brown PetscFunctionReturn(0); 939f1cd61daSJed Brown } 940f1cd61daSJed Brown 941f1cd61daSJed Brown #undef __FUNCT__ 9422bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps" 9432bfc04deSHong Zhang /*@ 9442bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 9452bfc04deSHong Zhang 9462bfc04deSHong Zhang Input Parameter: 9472bfc04deSHong Zhang + ts - the time-step context 9482bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 9492bfc04deSHong Zhang 9502bfc04deSHong Zhang Level: beginner 9512bfc04deSHong Zhang 9522bfc04deSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 9532bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 9542bfc04deSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 9552bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 9562bfc04deSHong Zhang @*/ 9577087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft) 9582bfc04deSHong Zhang { 9594ac538c5SBarry Smith PetscErrorCode ierr; 9602bfc04deSHong Zhang 9612bfc04deSHong Zhang PetscFunctionBegin; 9624ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 9632bfc04deSHong Zhang PetscFunctionReturn(0); 9642bfc04deSHong Zhang } 9657cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 9667cb94ee6SHong Zhang /*MC 96796f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 9687cb94ee6SHong Zhang 9697cb94ee6SHong Zhang Options Database: 9707cb94ee6SHong Zhang + -ts_sundials_type <bdf,adams> 9717cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 9727cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 9737cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 9747cb94ee6SHong Zhang . -ts_sundials_linear_tolerance <tol> 9757cb94ee6SHong Zhang . -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions 97616016685SHong Zhang . -ts_sundials_exact_final_time - Interpolate output to stop exactly at the final time 97716016685SHong Zhang - -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps 97816016685SHong Zhang 9797cb94ee6SHong Zhang 9807cb94ee6SHong Zhang Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply 9817cb94ee6SHong Zhang only PETSc PC options 9827cb94ee6SHong Zhang 9837cb94ee6SHong Zhang Level: beginner 9847cb94ee6SHong Zhang 9857cb94ee6SHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(), 9867cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime() 9877cb94ee6SHong Zhang 9887cb94ee6SHong Zhang M*/ 9897cb94ee6SHong Zhang EXTERN_C_BEGIN 9907cb94ee6SHong Zhang #undef __FUNCT__ 9917cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials" 9927087cfbeSBarry Smith PetscErrorCode TSCreate_Sundials(TS ts) 9937cb94ee6SHong Zhang { 9947cb94ee6SHong Zhang TS_Sundials *cvode; 9957cb94ee6SHong Zhang PetscErrorCode ierr; 9967cb94ee6SHong Zhang 9977cb94ee6SHong Zhang PetscFunctionBegin; 99817186662SBarry Smith if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for nonlinear problems"); 999*277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Sundials; 100028adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 100128adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 10027cb94ee6SHong Zhang ts->ops->setup = TSSetUp_Sundials_Nonlinear; 10037cb94ee6SHong Zhang ts->ops->step = TSStep_Sundials_Nonlinear; 10047cb94ee6SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_Sundials_Nonlinear; 10057cb94ee6SHong Zhang 100638f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr); 10077adad957SLisandro Dalcin ierr = PCCreate(((PetscObject)ts)->comm,&cvode->pc);CHKERRQ(ierr); 10087cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->pc);CHKERRQ(ierr); 10097cb94ee6SHong Zhang ts->data = (void*)cvode; 10106fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 10116fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 10127cb94ee6SHong Zhang cvode->restart = 5; 10137cb94ee6SHong Zhang cvode->linear_tol = .05; 10147cb94ee6SHong Zhang 10152bfc04deSHong Zhang cvode->exact_final_time = PETSC_TRUE; 10162bfc04deSHong Zhang cvode->monitorstep = PETSC_FALSE; 10177cb94ee6SHong Zhang 1018a07356b8SSatish Balay ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr); 1019f1cd61daSJed Brown 1020f1cd61daSJed Brown cvode->mindt = -1.; 1021f1cd61daSJed Brown cvode->maxdt = -1.; 1022f1cd61daSJed Brown 10237cb94ee6SHong Zhang /* set tolerance for Sundials */ 10247cb94ee6SHong Zhang cvode->reltol = 1e-6; 10252c823083SHong Zhang cvode->abstol = 1e-6; 10267cb94ee6SHong Zhang 10277cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials", 10287cb94ee6SHong Zhang TSSundialsSetType_Sundials);CHKERRQ(ierr); 10297cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C", 10307cb94ee6SHong Zhang "TSSundialsSetGMRESRestart_Sundials", 10317cb94ee6SHong Zhang TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr); 10327cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C", 10337cb94ee6SHong Zhang "TSSundialsSetLinearTolerance_Sundials", 10347cb94ee6SHong Zhang TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 10357cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C", 10367cb94ee6SHong Zhang "TSSundialsSetGramSchmidtType_Sundials", 10377cb94ee6SHong Zhang TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 10387cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C", 10397cb94ee6SHong Zhang "TSSundialsSetTolerance_Sundials", 10407cb94ee6SHong Zhang TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 1041f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C", 1042f1cd61daSJed Brown "TSSundialsSetMinTimeStep_Sundials", 1043f1cd61daSJed Brown TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr); 1044f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C", 1045f1cd61daSJed Brown "TSSundialsSetMaxTimeStep_Sundials", 1046f1cd61daSJed Brown TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr); 10477cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C", 10487cb94ee6SHong Zhang "TSSundialsGetPC_Sundials", 10497cb94ee6SHong Zhang TSSundialsGetPC_Sundials);CHKERRQ(ierr); 10507cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C", 10517cb94ee6SHong Zhang "TSSundialsGetIterations_Sundials", 10527cb94ee6SHong Zhang TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 10537cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C", 10547cb94ee6SHong Zhang "TSSundialsSetExactFinalTime_Sundials", 10557cb94ee6SHong Zhang TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr); 10562bfc04deSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C", 10572bfc04deSHong Zhang "TSSundialsMonitorInternalSteps_Sundials", 10582bfc04deSHong Zhang TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr); 10592bfc04deSHong Zhang 10607cb94ee6SHong Zhang PetscFunctionReturn(0); 10617cb94ee6SHong Zhang } 10627cb94ee6SHong Zhang EXTERN_C_END 1063