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; 1379f94935aSHong Zhang long int its,nsteps; 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 } 1569f94935aSHong Zhang 1579f94935aSHong Zhang if (flag){ /* display error message */ 1589f94935aSHong Zhang switch (flag){ 1599f94935aSHong Zhang case CV_ILL_INPUT: 1609f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT"); 1619f94935aSHong Zhang break; 1629f94935aSHong Zhang case CV_TOO_CLOSE: 1639f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE"); 1649f94935aSHong Zhang break; 1653c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1669f94935aSHong Zhang PetscReal tcur; 1679f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 1689f94935aSHong Zhang ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr); 1699f94935aSHong 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); 1703c7fefeeSJed Brown } break; 1719f94935aSHong Zhang case CV_TOO_MUCH_ACC: 1729f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC"); 1739f94935aSHong Zhang break; 1749f94935aSHong Zhang case CV_ERR_FAILURE: 1759f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE"); 1769f94935aSHong Zhang break; 1779f94935aSHong Zhang case CV_CONV_FAILURE: 1789f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE"); 1799f94935aSHong Zhang break; 1809f94935aSHong Zhang case CV_LINIT_FAIL: 1819f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL"); 1829f94935aSHong Zhang break; 1839f94935aSHong Zhang case CV_LSETUP_FAIL: 1849f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL"); 1859f94935aSHong Zhang break; 1869f94935aSHong Zhang case CV_LSOLVE_FAIL: 1879f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL"); 1889f94935aSHong Zhang break; 1899f94935aSHong Zhang case CV_RHSFUNC_FAIL: 1909f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL"); 1919f94935aSHong Zhang break; 1929f94935aSHong Zhang case CV_FIRST_RHSFUNC_ERR: 1939f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR"); 1949f94935aSHong Zhang break; 1959f94935aSHong Zhang case CV_REPTD_RHSFUNC_ERR: 1969f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR"); 1979f94935aSHong Zhang break; 1989f94935aSHong Zhang case CV_UNREC_RHSFUNC_ERR: 1999f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR"); 2009f94935aSHong Zhang break; 2019f94935aSHong Zhang case CV_RTFUNC_FAIL: 2029f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL"); 2039f94935aSHong Zhang break; 2049f94935aSHong Zhang default: 2059f94935aSHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag); 2069f94935aSHong Zhang } 2079f94935aSHong Zhang } 2089f94935aSHong Zhang 2097cb94ee6SHong Zhang if (t > ts->max_time && cvode->exact_final_time) { 2107cb94ee6SHong Zhang /* interpolate to final requested time */ 2117cb94ee6SHong Zhang ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr); 2127cb94ee6SHong Zhang t = tout; 2137cb94ee6SHong Zhang } 2147cb94ee6SHong Zhang ts->time_step = t - ts->ptime; 2157cb94ee6SHong Zhang ts->ptime = t; 2167cb94ee6SHong Zhang 2177cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 2187cb94ee6SHong Zhang ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr); 2197cb94ee6SHong Zhang ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr); 2207cb94ee6SHong Zhang ierr = VecResetArray(cvode->w1); CHKERRQ(ierr); 2217cb94ee6SHong Zhang ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr); 2227cb94ee6SHong Zhang ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr); 2237cb94ee6SHong Zhang ts->nonlinear_its = its; 2244ac7836bSHong Zhang ierr = CVSpilsGetNumLinIters(mem, &its); 2257cb94ee6SHong Zhang ts->linear_its = its; 2267cb94ee6SHong Zhang ts->steps++; 2273f2090d5SJed Brown ierr = TSPostStep(ts);CHKERRQ(ierr); 2287cb94ee6SHong Zhang ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr); 2297cb94ee6SHong Zhang } 2309f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 2319f94935aSHong Zhang *steps = nsteps; 2327cb94ee6SHong Zhang *time = t; 2337cb94ee6SHong Zhang PetscFunctionReturn(0); 2347cb94ee6SHong Zhang } 2357cb94ee6SHong Zhang 2367cb94ee6SHong Zhang #undef __FUNCT__ 2377cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials" 2387cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts) 2397cb94ee6SHong Zhang { 2407cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2417cb94ee6SHong Zhang PetscErrorCode ierr; 2427cb94ee6SHong Zhang 2437cb94ee6SHong Zhang PetscFunctionBegin; 2447cb94ee6SHong Zhang if (cvode->pmat) {ierr = MatDestroy(cvode->pmat);CHKERRQ(ierr);} 2457cb94ee6SHong Zhang if (cvode->pc) {ierr = PCDestroy(cvode->pc);CHKERRQ(ierr);} 2467cb94ee6SHong Zhang if (cvode->update) {ierr = VecDestroy(cvode->update);CHKERRQ(ierr);} 2477cb94ee6SHong Zhang if (cvode->func) {ierr = VecDestroy(cvode->func);CHKERRQ(ierr);} 2487cb94ee6SHong Zhang if (cvode->rhs) {ierr = VecDestroy(cvode->rhs);CHKERRQ(ierr);} 2497cb94ee6SHong Zhang if (cvode->w1) {ierr = VecDestroy(cvode->w1);CHKERRQ(ierr);} 2507cb94ee6SHong Zhang if (cvode->w2) {ierr = VecDestroy(cvode->w2);CHKERRQ(ierr);} 251a07356b8SSatish Balay ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr); 2522c823083SHong Zhang if (cvode->mem) {CVodeFree(&cvode->mem);} 2537cb94ee6SHong Zhang ierr = PetscFree(cvode);CHKERRQ(ierr); 2547cb94ee6SHong Zhang PetscFunctionReturn(0); 2557cb94ee6SHong Zhang } 2567cb94ee6SHong Zhang 2577cb94ee6SHong Zhang #undef __FUNCT__ 2587cb94ee6SHong Zhang #define __FUNCT__ "TSSetUp_Sundials_Nonlinear" 2597cb94ee6SHong Zhang PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts) 2607cb94ee6SHong Zhang { 2617cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2627cb94ee6SHong Zhang PetscErrorCode ierr; 26316016685SHong Zhang PetscInt glosize,locsize,i,flag; 2647cb94ee6SHong Zhang PetscScalar *y_data,*parray; 26516016685SHong Zhang void *mem; 26616016685SHong Zhang const PCType pctype; 267ace3abfcSBarry Smith PetscBool pcnone; 26816016685SHong Zhang Vec sol = ts->vec_sol; 2697cb94ee6SHong Zhang 2707cb94ee6SHong Zhang PetscFunctionBegin; 2717cb94ee6SHong Zhang ierr = PCSetFromOptions(cvode->pc);CHKERRQ(ierr); 2729f94935aSHong Zhang 2737cb94ee6SHong Zhang /* get the vector size */ 2747cb94ee6SHong Zhang ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr); 2757cb94ee6SHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr); 2767cb94ee6SHong Zhang 2777cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 278a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 279e32f2f54SBarry Smith if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated"); 2807cb94ee6SHong Zhang 28128adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 2827cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr); 2837cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y); 2847cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 2857cb94ee6SHong Zhang /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/ 2867cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 2877cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr); 2887cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->func);CHKERRQ(ierr); 2897cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr); 2907cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->func);CHKERRQ(ierr); 2917cb94ee6SHong Zhang 2927cb94ee6SHong Zhang /* 2937cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 2947cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 2957cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 2967cb94ee6SHong Zhang */ 2977adad957SLisandro Dalcin ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr); 2987adad957SLisandro Dalcin ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr); 2997cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr); 3007cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr); 30116016685SHong Zhang 30216016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 30316016685SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 304e32f2f54SBarry Smith if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails"); 30516016685SHong Zhang cvode->mem = mem; 30616016685SHong Zhang 30716016685SHong Zhang /* Set the pointer to user-defined data */ 30816016685SHong Zhang flag = CVodeSetUserData(mem, ts); 309e32f2f54SBarry Smith if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails"); 31016016685SHong Zhang 311fc6b9e64SJed Brown /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 312fc6b9e64SJed Brown flag = CVodeSetInitStep(mem,(realtype)ts->initial_time_step); 313fc6b9e64SJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed"); 314f1cd61daSJed Brown if (cvode->mindt > 0) { 315f1cd61daSJed Brown flag = CVodeSetMinStep(mem,(realtype)cvode->mindt); 3169f94935aSHong Zhang if (flag){ 3179f94935aSHong Zhang if (flag == CV_MEM_NULL){ 3189f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 3199f94935aSHong Zhang } else if (flag == CV_ILL_INPUT){ 3209f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size"); 3219f94935aSHong Zhang } else { 3229f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed"); 3239f94935aSHong Zhang } 3249f94935aSHong Zhang } 325f1cd61daSJed Brown } 326f1cd61daSJed Brown if (cvode->maxdt > 0) { 327f1cd61daSJed Brown flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt); 328f1cd61daSJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed"); 329f1cd61daSJed Brown } 330f1cd61daSJed Brown 33116016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 33216016685SHong Zhang * user's right hand side function in u'=f(t,u), the inital time T0, and 33316016685SHong Zhang * the initial dependent variable vector cvode->y */ 33416016685SHong Zhang flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 33516016685SHong Zhang if (flag){ 336e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag); 33716016685SHong Zhang } 33816016685SHong Zhang 3399f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 34016016685SHong Zhang flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 34116016685SHong Zhang if (flag){ 342e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag); 34316016685SHong Zhang } 34416016685SHong Zhang 3459f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 3469f94935aSHong Zhang flag = CVodeSetMaxNumSteps(mem,ts->max_steps); 34716016685SHong Zhang ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 34816016685SHong Zhang 34916016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 35016016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 35116016685SHong Zhang ierr = PCGetType(cvode->pc,&pctype);CHKERRQ(ierr); 35216016685SHong Zhang ierr = PetscTypeCompare((PetscObject)cvode->pc,PCNONE,&pcnone);CHKERRQ(ierr); 35316016685SHong Zhang if (pcnone){ 35416016685SHong Zhang flag = CVSpgmr(mem,PREC_NONE,0); 355e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 35616016685SHong Zhang } else { 35716016685SHong Zhang flag = CVSpgmr(mem,PREC_LEFT,0); 358e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 35916016685SHong Zhang 36016016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 36116016685SHong Zhang and the pointer to the user-defined block data */ 36216016685SHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 363e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 36416016685SHong Zhang } 36516016685SHong Zhang 36616016685SHong Zhang flag = CVSpilsSetGSType(mem, MODIFIED_GS); 367e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag); 3687cb94ee6SHong Zhang PetscFunctionReturn(0); 3697cb94ee6SHong Zhang } 3707cb94ee6SHong Zhang 3716fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 372dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0}; 3736fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 374dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0}; 375a04cf4d8SBarry Smith 3767cb94ee6SHong Zhang #undef __FUNCT__ 3777cb94ee6SHong Zhang #define __FUNCT__ "TSSetFromOptions_Sundials_Nonlinear" 3787cb94ee6SHong Zhang PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts) 3797cb94ee6SHong Zhang { 3807cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 3817cb94ee6SHong Zhang PetscErrorCode ierr; 3827cb94ee6SHong Zhang int indx; 383ace3abfcSBarry Smith PetscBool flag; 3847cb94ee6SHong Zhang 3857cb94ee6SHong Zhang PetscFunctionBegin; 3867cb94ee6SHong Zhang ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr); 3876fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr); 3887cb94ee6SHong Zhang if (flag) { 3896fadb2cdSHong Zhang ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr); 3907cb94ee6SHong Zhang } 3916fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr); 3927cb94ee6SHong Zhang if (flag) { 3937cb94ee6SHong Zhang ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 3947cb94ee6SHong Zhang } 3957cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr); 3967cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr); 397f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr); 398f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr); 3997cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr); 4007cb94ee6SHong Zhang ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr); 401*acfcf0e5SJed 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); 402*acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr); 4037cb94ee6SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4047cb94ee6SHong Zhang PetscFunctionReturn(0); 4057cb94ee6SHong Zhang } 4067cb94ee6SHong Zhang 4077cb94ee6SHong Zhang #undef __FUNCT__ 4087cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials" 4097cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 4107cb94ee6SHong Zhang { 4117cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4127cb94ee6SHong Zhang PetscErrorCode ierr; 4137cb94ee6SHong Zhang char *type; 4147cb94ee6SHong Zhang char atype[] = "Adams"; 4157cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 416ace3abfcSBarry Smith PetscBool iascii,isstring; 4172c823083SHong Zhang long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 4182c823083SHong Zhang PetscInt qlast,qcur; 4195d47aee6SHong Zhang PetscReal hinused,hlast,hcur,tcur,tolsfac; 4207cb94ee6SHong Zhang 4217cb94ee6SHong Zhang PetscFunctionBegin; 4227cb94ee6SHong Zhang if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;} 4237cb94ee6SHong Zhang else {type = btype;} 4247cb94ee6SHong Zhang 4252692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4262692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 4277cb94ee6SHong Zhang if (iascii) { 4287cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr); 4297cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr); 4307cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr); 4317cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr); 4327cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr); 4337cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 4347cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4357cb94ee6SHong Zhang } else { 4367cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4377cb94ee6SHong Zhang } 438f1cd61daSJed Brown if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);} 439f1cd61daSJed Brown if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);} 4402c823083SHong Zhang 4415d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4425d47aee6SHong Zhang ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr); 4435d47aee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr); 4442c823083SHong Zhang ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals, 4452c823083SHong Zhang &nlinsetups,&nfails,&qlast,&qcur, 4462c823083SHong Zhang &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr); 4472c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr); 4482c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr); 4492c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr); 4502c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr); 4512c823083SHong Zhang 4522c823083SHong Zhang ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr); 4532c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr); 4542c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr); 4552c823083SHong Zhang 4562c823083SHong Zhang ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */ 4572c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr); 4582c823083SHong Zhang ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr); 4592c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr); 4602c823083SHong Zhang ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4612c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr); 4622c823083SHong Zhang ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr); 4632c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr); 4642c823083SHong Zhang ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4652c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr); 4662c823083SHong Zhang ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4672c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr); 4687cb94ee6SHong Zhang } else if (isstring) { 4697cb94ee6SHong Zhang ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 4707cb94ee6SHong Zhang } else { 471e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name); 4727cb94ee6SHong Zhang } 4737cb94ee6SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 4747cb94ee6SHong Zhang ierr = PCView(cvode->pc,viewer);CHKERRQ(ierr); 4757cb94ee6SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 4767cb94ee6SHong Zhang PetscFunctionReturn(0); 4777cb94ee6SHong Zhang } 4787cb94ee6SHong Zhang 4797cb94ee6SHong Zhang 4807cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 4817cb94ee6SHong Zhang EXTERN_C_BEGIN 4827cb94ee6SHong Zhang #undef __FUNCT__ 4837cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials" 4846fadb2cdSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 4857cb94ee6SHong Zhang { 4867cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4877cb94ee6SHong Zhang 4887cb94ee6SHong Zhang PetscFunctionBegin; 4897cb94ee6SHong Zhang cvode->cvode_type = type; 4907cb94ee6SHong Zhang PetscFunctionReturn(0); 4917cb94ee6SHong Zhang } 4927cb94ee6SHong Zhang EXTERN_C_END 4937cb94ee6SHong Zhang 4947cb94ee6SHong Zhang EXTERN_C_BEGIN 4957cb94ee6SHong Zhang #undef __FUNCT__ 4967cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials" 4977cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart_Sundials(TS ts,int restart) 4987cb94ee6SHong Zhang { 4997cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5007cb94ee6SHong Zhang 5017cb94ee6SHong Zhang PetscFunctionBegin; 5027cb94ee6SHong Zhang cvode->restart = restart; 5037cb94ee6SHong Zhang PetscFunctionReturn(0); 5047cb94ee6SHong Zhang } 5057cb94ee6SHong Zhang EXTERN_C_END 5067cb94ee6SHong Zhang 5077cb94ee6SHong Zhang EXTERN_C_BEGIN 5087cb94ee6SHong Zhang #undef __FUNCT__ 5097cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials" 5107cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance_Sundials(TS ts,double tol) 5117cb94ee6SHong Zhang { 5127cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5137cb94ee6SHong Zhang 5147cb94ee6SHong Zhang PetscFunctionBegin; 5157cb94ee6SHong Zhang cvode->linear_tol = tol; 5167cb94ee6SHong Zhang PetscFunctionReturn(0); 5177cb94ee6SHong Zhang } 5187cb94ee6SHong Zhang EXTERN_C_END 5197cb94ee6SHong Zhang 5207cb94ee6SHong Zhang EXTERN_C_BEGIN 5217cb94ee6SHong Zhang #undef __FUNCT__ 5227cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials" 5237cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 5247cb94ee6SHong Zhang { 5257cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5267cb94ee6SHong Zhang 5277cb94ee6SHong Zhang PetscFunctionBegin; 5287cb94ee6SHong Zhang cvode->gtype = type; 5297cb94ee6SHong Zhang PetscFunctionReturn(0); 5307cb94ee6SHong Zhang } 5317cb94ee6SHong Zhang EXTERN_C_END 5327cb94ee6SHong Zhang 5337cb94ee6SHong Zhang EXTERN_C_BEGIN 5347cb94ee6SHong Zhang #undef __FUNCT__ 5357cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials" 5367cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel) 5377cb94ee6SHong Zhang { 5387cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5397cb94ee6SHong Zhang 5407cb94ee6SHong Zhang PetscFunctionBegin; 5417cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 5427cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 5437cb94ee6SHong Zhang PetscFunctionReturn(0); 5447cb94ee6SHong Zhang } 5457cb94ee6SHong Zhang EXTERN_C_END 5467cb94ee6SHong Zhang 5477cb94ee6SHong Zhang EXTERN_C_BEGIN 5487cb94ee6SHong Zhang #undef __FUNCT__ 549f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials" 550f1cd61daSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 551f1cd61daSJed Brown { 552f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 553f1cd61daSJed Brown 554f1cd61daSJed Brown PetscFunctionBegin; 555f1cd61daSJed Brown cvode->mindt = mindt; 556f1cd61daSJed Brown PetscFunctionReturn(0); 557f1cd61daSJed Brown } 558f1cd61daSJed Brown EXTERN_C_END 559f1cd61daSJed Brown 560f1cd61daSJed Brown EXTERN_C_BEGIN 561f1cd61daSJed Brown #undef __FUNCT__ 562f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials" 563f1cd61daSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 564f1cd61daSJed Brown { 565f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 566f1cd61daSJed Brown 567f1cd61daSJed Brown PetscFunctionBegin; 568f1cd61daSJed Brown cvode->maxdt = maxdt; 569f1cd61daSJed Brown PetscFunctionReturn(0); 570f1cd61daSJed Brown } 571f1cd61daSJed Brown EXTERN_C_END 572f1cd61daSJed Brown 573f1cd61daSJed Brown EXTERN_C_BEGIN 574f1cd61daSJed Brown #undef __FUNCT__ 5757cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials" 5767cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC_Sundials(TS ts,PC *pc) 5777cb94ee6SHong Zhang { 5787cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5797cb94ee6SHong Zhang 5807cb94ee6SHong Zhang PetscFunctionBegin; 5817cb94ee6SHong Zhang *pc = cvode->pc; 5827cb94ee6SHong Zhang PetscFunctionReturn(0); 5837cb94ee6SHong Zhang } 5847cb94ee6SHong Zhang EXTERN_C_END 5857cb94ee6SHong Zhang 5867cb94ee6SHong Zhang EXTERN_C_BEGIN 5877cb94ee6SHong Zhang #undef __FUNCT__ 5887cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials" 5897cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 5907cb94ee6SHong Zhang { 5917cb94ee6SHong Zhang PetscFunctionBegin; 5922c823083SHong Zhang if (nonlin) *nonlin = ts->nonlinear_its; 5932c823083SHong Zhang if (lin) *lin = ts->linear_its; 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__ "TSSundialsSetExactFinalTime_Sundials" 601ace3abfcSBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime_Sundials(TS ts,PetscBool s) 6027cb94ee6SHong Zhang { 6037cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 6047cb94ee6SHong Zhang 6057cb94ee6SHong Zhang PetscFunctionBegin; 6067cb94ee6SHong Zhang cvode->exact_final_time = s; 6077cb94ee6SHong Zhang PetscFunctionReturn(0); 6087cb94ee6SHong Zhang } 6097cb94ee6SHong Zhang EXTERN_C_END 6102bfc04deSHong Zhang 6112bfc04deSHong Zhang EXTERN_C_BEGIN 6122bfc04deSHong Zhang #undef __FUNCT__ 6132bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials" 614ace3abfcSBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s) 6152bfc04deSHong Zhang { 6162bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 6172bfc04deSHong Zhang 6182bfc04deSHong Zhang PetscFunctionBegin; 6192bfc04deSHong Zhang cvode->monitorstep = s; 6202bfc04deSHong Zhang PetscFunctionReturn(0); 6212bfc04deSHong Zhang } 6222bfc04deSHong Zhang EXTERN_C_END 6237cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 6247cb94ee6SHong Zhang 6257cb94ee6SHong Zhang #undef __FUNCT__ 6267cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations" 6277cb94ee6SHong Zhang /*@C 6287cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 6297cb94ee6SHong Zhang 6307cb94ee6SHong Zhang Not Collective 6317cb94ee6SHong Zhang 6327cb94ee6SHong Zhang Input parameters: 6337cb94ee6SHong Zhang . ts - the time-step context 6347cb94ee6SHong Zhang 6357cb94ee6SHong Zhang Output Parameters: 6367cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 6377cb94ee6SHong Zhang - lin - number of linear iterations 6387cb94ee6SHong Zhang 6397cb94ee6SHong Zhang Level: advanced 6407cb94ee6SHong Zhang 6417cb94ee6SHong Zhang Notes: 6427cb94ee6SHong Zhang These return the number since the creation of the TS object 6437cb94ee6SHong Zhang 6447cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations 6457cb94ee6SHong Zhang 6467cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6477cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 6487cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 649f1cd61daSJed Brown TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSundialsSetExactFinalTime() 6507cb94ee6SHong Zhang 6517cb94ee6SHong Zhang @*/ 6527cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 6537cb94ee6SHong Zhang { 6544ac538c5SBarry Smith PetscErrorCode ierr; 6557cb94ee6SHong Zhang 6567cb94ee6SHong Zhang PetscFunctionBegin; 6574ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr); 6587cb94ee6SHong Zhang PetscFunctionReturn(0); 6597cb94ee6SHong Zhang } 6607cb94ee6SHong Zhang 6617cb94ee6SHong Zhang #undef __FUNCT__ 6627cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType" 6637cb94ee6SHong Zhang /*@ 6647cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 6657cb94ee6SHong Zhang 6663f9fe445SBarry Smith Logically Collective on TS 6677cb94ee6SHong Zhang 6687cb94ee6SHong Zhang Input parameters: 6697cb94ee6SHong Zhang + ts - the time-step context 6707cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 6717cb94ee6SHong Zhang 6727cb94ee6SHong Zhang Level: intermediate 6737cb94ee6SHong Zhang 6747cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula 6757cb94ee6SHong Zhang 6767cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetGMRESRestart(), 6777cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 6787cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6797cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 6807cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 6817cb94ee6SHong Zhang @*/ 6826fadb2cdSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType(TS ts,TSSundialsLmmType type) 6837cb94ee6SHong Zhang { 6844ac538c5SBarry Smith PetscErrorCode ierr; 6857cb94ee6SHong Zhang 6867cb94ee6SHong Zhang PetscFunctionBegin; 6874ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr); 6887cb94ee6SHong Zhang PetscFunctionReturn(0); 6897cb94ee6SHong Zhang } 6907cb94ee6SHong Zhang 6917cb94ee6SHong Zhang #undef __FUNCT__ 6927cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart" 6937cb94ee6SHong Zhang /*@ 6947cb94ee6SHong Zhang TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by 6957cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 6967cb94ee6SHong Zhang this is ALSO the maximum number of GMRES steps that will be used. 6977cb94ee6SHong Zhang 6983f9fe445SBarry Smith Logically Collective on TS 6997cb94ee6SHong Zhang 7007cb94ee6SHong Zhang Input parameters: 7017cb94ee6SHong Zhang + ts - the time-step context 7027cb94ee6SHong Zhang - restart - number of direction vectors (the restart size). 7037cb94ee6SHong Zhang 7047cb94ee6SHong Zhang Level: advanced 7057cb94ee6SHong Zhang 7067cb94ee6SHong Zhang .keywords: GMRES, restart 7077cb94ee6SHong Zhang 7087cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 7097cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 7107cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7117cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 7127cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 7137cb94ee6SHong Zhang 7147cb94ee6SHong Zhang @*/ 7157cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart(TS ts,int restart) 7167cb94ee6SHong Zhang { 7174ac538c5SBarry Smith PetscErrorCode ierr; 7187cb94ee6SHong Zhang 7197cb94ee6SHong Zhang PetscFunctionBegin; 720c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(ts,restart,2); 7214ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetGMRESRestart_C",(TS,int),(ts,restart));CHKERRQ(ierr); 7227cb94ee6SHong Zhang PetscFunctionReturn(0); 7237cb94ee6SHong Zhang } 7247cb94ee6SHong Zhang 7257cb94ee6SHong Zhang #undef __FUNCT__ 7267cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance" 7277cb94ee6SHong Zhang /*@ 7287cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 7297cb94ee6SHong Zhang system by SUNDIALS. 7307cb94ee6SHong Zhang 7313f9fe445SBarry Smith Logically Collective on TS 7327cb94ee6SHong Zhang 7337cb94ee6SHong Zhang Input parameters: 7347cb94ee6SHong Zhang + ts - the time-step context 7357cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 7367cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 7377cb94ee6SHong Zhang 7387cb94ee6SHong Zhang Level: advanced 7397cb94ee6SHong Zhang 7407cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS 7417cb94ee6SHong Zhang 7427cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7437cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 7447cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7457cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 7467cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 7477cb94ee6SHong Zhang 7487cb94ee6SHong Zhang @*/ 7497cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance(TS ts,double tol) 7507cb94ee6SHong Zhang { 7514ac538c5SBarry Smith PetscErrorCode ierr; 7527cb94ee6SHong Zhang 7537cb94ee6SHong Zhang PetscFunctionBegin; 754c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,tol,2); 7554ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr); 7567cb94ee6SHong Zhang PetscFunctionReturn(0); 7577cb94ee6SHong Zhang } 7587cb94ee6SHong Zhang 7597cb94ee6SHong Zhang #undef __FUNCT__ 7607cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType" 7617cb94ee6SHong Zhang /*@ 7627cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 7637cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 7647cb94ee6SHong Zhang 7653f9fe445SBarry Smith Logically Collective on TS 7667cb94ee6SHong Zhang 7677cb94ee6SHong Zhang Input parameters: 7687cb94ee6SHong Zhang + ts - the time-step context 7697cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 7707cb94ee6SHong Zhang 7717cb94ee6SHong Zhang Level: advanced 7727cb94ee6SHong Zhang 7737cb94ee6SHong Zhang .keywords: Sundials, orthogonalization 7747cb94ee6SHong Zhang 7757cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7767cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 7777cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7787cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 7797cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 7807cb94ee6SHong Zhang 7817cb94ee6SHong Zhang @*/ 7827cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 7837cb94ee6SHong Zhang { 7844ac538c5SBarry Smith PetscErrorCode ierr; 7857cb94ee6SHong Zhang 7867cb94ee6SHong Zhang PetscFunctionBegin; 7874ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr); 7887cb94ee6SHong Zhang PetscFunctionReturn(0); 7897cb94ee6SHong Zhang } 7907cb94ee6SHong Zhang 7917cb94ee6SHong Zhang #undef __FUNCT__ 7927cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance" 7937cb94ee6SHong Zhang /*@ 7947cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 7957cb94ee6SHong Zhang Sundials for error control. 7967cb94ee6SHong Zhang 7973f9fe445SBarry Smith Logically Collective on TS 7987cb94ee6SHong Zhang 7997cb94ee6SHong Zhang Input parameters: 8007cb94ee6SHong Zhang + ts - the time-step context 8017cb94ee6SHong Zhang . aabs - the absolute tolerance 8027cb94ee6SHong Zhang - rel - the relative tolerance 8037cb94ee6SHong Zhang 8047cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 8057cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 8067cb94ee6SHong Zhang 8077cb94ee6SHong Zhang Level: intermediate 8087cb94ee6SHong Zhang 8097cb94ee6SHong Zhang .keywords: Sundials, tolerance 8107cb94ee6SHong Zhang 8117cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8127cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 8137cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8147cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 8157cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 8167cb94ee6SHong Zhang 8177cb94ee6SHong Zhang @*/ 8187cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance(TS ts,double aabs,double rel) 8197cb94ee6SHong Zhang { 8204ac538c5SBarry Smith PetscErrorCode ierr; 8217cb94ee6SHong Zhang 8227cb94ee6SHong Zhang PetscFunctionBegin; 8234ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr); 8247cb94ee6SHong Zhang PetscFunctionReturn(0); 8257cb94ee6SHong Zhang } 8267cb94ee6SHong Zhang 8277cb94ee6SHong Zhang #undef __FUNCT__ 8287cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC" 8297cb94ee6SHong Zhang /*@ 8307cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 8317cb94ee6SHong Zhang 8327cb94ee6SHong Zhang Input Parameter: 8337cb94ee6SHong Zhang . ts - the time-step context 8347cb94ee6SHong Zhang 8357cb94ee6SHong Zhang Output Parameter: 8367cb94ee6SHong Zhang . pc - the preconditioner context 8377cb94ee6SHong Zhang 8387cb94ee6SHong Zhang Level: advanced 8397cb94ee6SHong Zhang 8407cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8417cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 8427cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8437cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 8447cb94ee6SHong Zhang @*/ 8457cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC(TS ts,PC *pc) 8467cb94ee6SHong Zhang { 8474ac538c5SBarry Smith PetscErrorCode ierr; 8487cb94ee6SHong Zhang 8497cb94ee6SHong Zhang PetscFunctionBegin; 8504ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr); 8517cb94ee6SHong Zhang PetscFunctionReturn(0); 8527cb94ee6SHong Zhang } 8537cb94ee6SHong Zhang 8547cb94ee6SHong Zhang #undef __FUNCT__ 8557cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime" 8567cb94ee6SHong Zhang /*@ 8577cb94ee6SHong Zhang TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the 8587cb94ee6SHong Zhang exact final time requested by the user or just returns it at the final time 8597cb94ee6SHong Zhang it computed. (Defaults to true). 8607cb94ee6SHong Zhang 8617cb94ee6SHong Zhang Input Parameter: 8627cb94ee6SHong Zhang + ts - the time-step context 8637cb94ee6SHong Zhang - ft - PETSC_TRUE if interpolates, else PETSC_FALSE 8647cb94ee6SHong Zhang 8657cb94ee6SHong Zhang Level: beginner 8667cb94ee6SHong Zhang 8677cb94ee6SHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8687cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 8697cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8707cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 8717cb94ee6SHong Zhang @*/ 872ace3abfcSBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime(TS ts,PetscBool ft) 8737cb94ee6SHong Zhang { 8744ac538c5SBarry Smith PetscErrorCode ierr; 8757cb94ee6SHong Zhang 8767cb94ee6SHong Zhang PetscFunctionBegin; 8774ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetExactFinalTime_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 8787cb94ee6SHong Zhang PetscFunctionReturn(0); 8797cb94ee6SHong Zhang } 8807cb94ee6SHong Zhang 8812bfc04deSHong Zhang #undef __FUNCT__ 882f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep" 883f1cd61daSJed Brown /*@ 884f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 885f1cd61daSJed Brown 886f1cd61daSJed Brown Input Parameter: 887f1cd61daSJed Brown + ts - the time-step context 888f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 889f1cd61daSJed Brown 890fc6b9e64SJed Brown Note: 891fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 892fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 893fc6b9e64SJed Brown 894f1cd61daSJed Brown Level: beginner 895f1cd61daSJed Brown 896f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 897f1cd61daSJed Brown @*/ 898f1cd61daSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 899f1cd61daSJed Brown { 9004ac538c5SBarry Smith PetscErrorCode ierr; 901f1cd61daSJed Brown 902f1cd61daSJed Brown PetscFunctionBegin; 9034ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr); 904f1cd61daSJed Brown PetscFunctionReturn(0); 905f1cd61daSJed Brown } 906f1cd61daSJed Brown 907f1cd61daSJed Brown #undef __FUNCT__ 908f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep" 909f1cd61daSJed Brown /*@ 910f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 911f1cd61daSJed Brown 912f1cd61daSJed Brown Input Parameter: 913f1cd61daSJed Brown + ts - the time-step context 914f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 915f1cd61daSJed Brown 916f1cd61daSJed Brown Level: beginner 917f1cd61daSJed Brown 918f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 919f1cd61daSJed Brown @*/ 920f1cd61daSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 921f1cd61daSJed Brown { 9224ac538c5SBarry Smith PetscErrorCode ierr; 923f1cd61daSJed Brown 924f1cd61daSJed Brown PetscFunctionBegin; 9254ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 926f1cd61daSJed Brown PetscFunctionReturn(0); 927f1cd61daSJed Brown } 928f1cd61daSJed Brown 929f1cd61daSJed Brown #undef __FUNCT__ 9302bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps" 9312bfc04deSHong Zhang /*@ 9322bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 9332bfc04deSHong Zhang 9342bfc04deSHong Zhang Input Parameter: 9352bfc04deSHong Zhang + ts - the time-step context 9362bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 9372bfc04deSHong Zhang 9382bfc04deSHong Zhang Level: beginner 9392bfc04deSHong Zhang 9402bfc04deSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 9412bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 9422bfc04deSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 9432bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 9442bfc04deSHong Zhang @*/ 945ace3abfcSBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSSundialsMonitorInternalSteps(TS ts,PetscBool ft) 9462bfc04deSHong Zhang { 9474ac538c5SBarry Smith PetscErrorCode ierr; 9482bfc04deSHong Zhang 9492bfc04deSHong Zhang PetscFunctionBegin; 9504ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 9512bfc04deSHong Zhang PetscFunctionReturn(0); 9522bfc04deSHong Zhang } 9537cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 9547cb94ee6SHong Zhang /*MC 95596f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 9567cb94ee6SHong Zhang 9577cb94ee6SHong Zhang Options Database: 9587cb94ee6SHong Zhang + -ts_sundials_type <bdf,adams> 9597cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 9607cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 9617cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 9627cb94ee6SHong Zhang . -ts_sundials_linear_tolerance <tol> 9637cb94ee6SHong Zhang . -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions 96416016685SHong Zhang . -ts_sundials_exact_final_time - Interpolate output to stop exactly at the final time 96516016685SHong Zhang - -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps 96616016685SHong Zhang 9677cb94ee6SHong Zhang 9687cb94ee6SHong Zhang Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply 9697cb94ee6SHong Zhang only PETSc PC options 9707cb94ee6SHong Zhang 9717cb94ee6SHong Zhang Level: beginner 9727cb94ee6SHong Zhang 9737cb94ee6SHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(), 9747cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime() 9757cb94ee6SHong Zhang 9767cb94ee6SHong Zhang M*/ 9777cb94ee6SHong Zhang EXTERN_C_BEGIN 9787cb94ee6SHong Zhang #undef __FUNCT__ 9797cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials" 9807cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Sundials(TS ts) 9817cb94ee6SHong Zhang { 9827cb94ee6SHong Zhang TS_Sundials *cvode; 9837cb94ee6SHong Zhang PetscErrorCode ierr; 9847cb94ee6SHong Zhang 9857cb94ee6SHong Zhang PetscFunctionBegin; 98617186662SBarry Smith if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for nonlinear problems"); 98728adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 98828adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 9897cb94ee6SHong Zhang ts->ops->setup = TSSetUp_Sundials_Nonlinear; 9907cb94ee6SHong Zhang ts->ops->step = TSStep_Sundials_Nonlinear; 9917cb94ee6SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_Sundials_Nonlinear; 9927cb94ee6SHong Zhang 99338f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr); 9947adad957SLisandro Dalcin ierr = PCCreate(((PetscObject)ts)->comm,&cvode->pc);CHKERRQ(ierr); 9957cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->pc);CHKERRQ(ierr); 9967cb94ee6SHong Zhang ts->data = (void*)cvode; 9976fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 9986fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 9997cb94ee6SHong Zhang cvode->restart = 5; 10007cb94ee6SHong Zhang cvode->linear_tol = .05; 10017cb94ee6SHong Zhang 10022bfc04deSHong Zhang cvode->exact_final_time = PETSC_TRUE; 10032bfc04deSHong Zhang cvode->monitorstep = PETSC_FALSE; 10047cb94ee6SHong Zhang 1005a07356b8SSatish Balay ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr); 1006f1cd61daSJed Brown 1007f1cd61daSJed Brown cvode->mindt = -1.; 1008f1cd61daSJed Brown cvode->maxdt = -1.; 1009f1cd61daSJed Brown 10107cb94ee6SHong Zhang /* set tolerance for Sundials */ 10117cb94ee6SHong Zhang cvode->reltol = 1e-6; 10122c823083SHong Zhang cvode->abstol = 1e-6; 10137cb94ee6SHong Zhang 10147cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials", 10157cb94ee6SHong Zhang TSSundialsSetType_Sundials);CHKERRQ(ierr); 10167cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C", 10177cb94ee6SHong Zhang "TSSundialsSetGMRESRestart_Sundials", 10187cb94ee6SHong Zhang TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr); 10197cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C", 10207cb94ee6SHong Zhang "TSSundialsSetLinearTolerance_Sundials", 10217cb94ee6SHong Zhang TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 10227cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C", 10237cb94ee6SHong Zhang "TSSundialsSetGramSchmidtType_Sundials", 10247cb94ee6SHong Zhang TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 10257cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C", 10267cb94ee6SHong Zhang "TSSundialsSetTolerance_Sundials", 10277cb94ee6SHong Zhang TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 1028f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C", 1029f1cd61daSJed Brown "TSSundialsSetMinTimeStep_Sundials", 1030f1cd61daSJed Brown TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr); 1031f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C", 1032f1cd61daSJed Brown "TSSundialsSetMaxTimeStep_Sundials", 1033f1cd61daSJed Brown TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr); 10347cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C", 10357cb94ee6SHong Zhang "TSSundialsGetPC_Sundials", 10367cb94ee6SHong Zhang TSSundialsGetPC_Sundials);CHKERRQ(ierr); 10377cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C", 10387cb94ee6SHong Zhang "TSSundialsGetIterations_Sundials", 10397cb94ee6SHong Zhang TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 10407cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C", 10417cb94ee6SHong Zhang "TSSundialsSetExactFinalTime_Sundials", 10427cb94ee6SHong Zhang TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr); 10432bfc04deSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C", 10442bfc04deSHong Zhang "TSSundialsMonitorInternalSteps_Sundials", 10452bfc04deSHong Zhang TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr); 10462bfc04deSHong Zhang 10477cb94ee6SHong Zhang PetscFunctionReturn(0); 10487cb94ee6SHong Zhang } 10497cb94ee6SHong Zhang EXTERN_C_END 1050