17cb94ee6SHong Zhang /* 2db268bfaSHong Zhang Provides a PETSc interface to SUNDIALS/CVODE solver. 37cb94ee6SHong Zhang The interface to PVODE (old version of CVODE) was originally contributed 47cb94ee6SHong Zhang by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik. 528adb3f7SHong Zhang 628adb3f7SHong Zhang Reference: sundials-2.4.0/examples/cvode/parallel/cvDiurnal_kry_p.c 77cb94ee6SHong Zhang */ 856a740aaSHong Zhang #include "sundials.h" /*I "petscts.h" I*/ 97cb94ee6SHong Zhang 107cb94ee6SHong Zhang /* 117cb94ee6SHong Zhang TSPrecond_Sundials - function that we provide to SUNDIALS to 127cb94ee6SHong Zhang evaluate the preconditioner. 137cb94ee6SHong Zhang */ 147cb94ee6SHong Zhang #undef __FUNCT__ 157cb94ee6SHong Zhang #define __FUNCT__ "TSPrecond_Sundials" 167cb94ee6SHong Zhang PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy, 177cb94ee6SHong Zhang booleantype jok,booleantype *jcurPtr, 187cb94ee6SHong Zhang realtype _gamma,void *P_data, 197cb94ee6SHong Zhang N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3) 207cb94ee6SHong Zhang { 217cb94ee6SHong Zhang TS ts = (TS) P_data; 227cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 23dcbc6d53SSean Farley PC pc; 247cb94ee6SHong Zhang PetscErrorCode ierr; 25*0679a0aeSJed Brown Mat J,P; 26*0679a0aeSJed Brown Vec yy = cvode->w1,yydot = cvode->ydot; 27*0679a0aeSJed Brown PetscReal gm = (PetscReal)_gamma; 287cb94ee6SHong Zhang MatStructure str = DIFFERENT_NONZERO_PATTERN; 297cb94ee6SHong Zhang PetscScalar *y_data; 307cb94ee6SHong Zhang 317cb94ee6SHong Zhang PetscFunctionBegin; 32*0679a0aeSJed Brown ierr = TSGetIJacobian(ts,&J,&P,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 337cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(y); 347cb94ee6SHong Zhang ierr = VecPlaceArray(yy,y_data); CHKERRQ(ierr); 35*0679a0aeSJed Brown ierr = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */ 36*0679a0aeSJed Brown /* compute the shifted Jacobian (1/gm)*I + Jrest */ 37*0679a0aeSJed Brown ierr = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,&J,&P,&str,PETSC_FALSE);CHKERRQ(ierr); 387cb94ee6SHong Zhang ierr = VecResetArray(yy); CHKERRQ(ierr); 39*0679a0aeSJed Brown ierr = MatScale(P,gm);CHKERRQ(ierr); /* turn into I-gm*Jrest, J is not used by Sundials */ 407cb94ee6SHong Zhang *jcurPtr = TRUE; 41dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr); 42*0679a0aeSJed Brown ierr = PCSetOperators(pc,J,P,str);CHKERRQ(ierr); 437cb94ee6SHong Zhang PetscFunctionReturn(0); 447cb94ee6SHong Zhang } 457cb94ee6SHong Zhang 467cb94ee6SHong Zhang /* 477cb94ee6SHong Zhang TSPSolve_Sundials - routine that we provide to Sundials that applies the preconditioner. 487cb94ee6SHong Zhang */ 497cb94ee6SHong Zhang #undef __FUNCT__ 507cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials" 514ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z, 524ac7836bSHong Zhang realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp) 537cb94ee6SHong Zhang { 547cb94ee6SHong Zhang TS ts = (TS) P_data; 557cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 56dcbc6d53SSean Farley PC pc; 577cb94ee6SHong Zhang Vec rr = cvode->w1,zz = cvode->w2; 587cb94ee6SHong Zhang PetscErrorCode ierr; 597cb94ee6SHong Zhang PetscScalar *r_data,*z_data; 607cb94ee6SHong Zhang 617cb94ee6SHong Zhang PetscFunctionBegin; 627cb94ee6SHong Zhang /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/ 637cb94ee6SHong Zhang r_data = (PetscScalar *) N_VGetArrayPointer(r); 647cb94ee6SHong Zhang z_data = (PetscScalar *) N_VGetArrayPointer(z); 657cb94ee6SHong Zhang ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr); 667cb94ee6SHong Zhang ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr); 674ac7836bSHong Zhang 687cb94ee6SHong Zhang /* Solve the Px=r and put the result in zz */ 69dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr); 707cb94ee6SHong Zhang ierr = PCApply(pc,rr,zz); CHKERRQ(ierr); 717cb94ee6SHong Zhang ierr = VecResetArray(rr); CHKERRQ(ierr); 727cb94ee6SHong Zhang ierr = VecResetArray(zz); CHKERRQ(ierr); 737cb94ee6SHong Zhang PetscFunctionReturn(0); 747cb94ee6SHong Zhang } 757cb94ee6SHong Zhang 767cb94ee6SHong Zhang /* 777cb94ee6SHong Zhang TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side. 787cb94ee6SHong Zhang */ 797cb94ee6SHong Zhang #undef __FUNCT__ 807cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials" 814ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx) 827cb94ee6SHong Zhang { 837cb94ee6SHong Zhang TS ts = (TS) ctx; 847adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)ts)->comm; 857cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 86*0679a0aeSJed Brown Vec yy = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot; 877cb94ee6SHong Zhang PetscScalar *y_data,*ydot_data; 887cb94ee6SHong Zhang PetscErrorCode ierr; 897cb94ee6SHong Zhang 907cb94ee6SHong Zhang PetscFunctionBegin; 915ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 927cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(y); 937cb94ee6SHong Zhang ydot_data = (PetscScalar *) N_VGetArrayPointer(ydot); 944ac7836bSHong Zhang ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr); 954ac7836bSHong Zhang ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr); 96*0679a0aeSJed Brown ierr = VecZeroEntries(yydot);CHKERRQ(ierr); 974ac7836bSHong Zhang 987cb94ee6SHong Zhang /* now compute the right hand side function */ 99*0679a0aeSJed Brown ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE); CHKERRABORT(comm,ierr); 100*0679a0aeSJed Brown ierr = VecScale(yyd,-1.);CHKERRQ(ierr); 1017cb94ee6SHong Zhang ierr = VecResetArray(yy); CHKERRABORT(comm,ierr); 1027cb94ee6SHong Zhang ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr); 1034ac7836bSHong Zhang PetscFunctionReturn(0); 1047cb94ee6SHong Zhang } 1057cb94ee6SHong Zhang 1067cb94ee6SHong Zhang /* 107214bc6a2SJed Brown TSStep_Sundials - Calls Sundials to integrate the ODE. 1087cb94ee6SHong Zhang */ 1097cb94ee6SHong Zhang #undef __FUNCT__ 110214bc6a2SJed Brown #define __FUNCT__ "TSStep_Sundials" 111193ac0bcSJed Brown PetscErrorCode TSStep_Sundials(TS ts) 1127cb94ee6SHong Zhang { 1137cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1147cb94ee6SHong Zhang Vec sol = ts->vec_sol; 1157cb94ee6SHong Zhang PetscErrorCode ierr; 116186e87acSLisandro Dalcin PetscInt i,flag; 1179f94935aSHong Zhang long int its,nsteps; 1187cb94ee6SHong Zhang realtype t,tout; 1197cb94ee6SHong Zhang PetscScalar *y_data; 1207cb94ee6SHong Zhang void *mem; 1217cb94ee6SHong Zhang 1227cb94ee6SHong Zhang PetscFunctionBegin; 12316016685SHong Zhang mem = cvode->mem; 1247cb94ee6SHong Zhang tout = ts->max_time; 1257cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr); 1267cb94ee6SHong Zhang N_VSetArrayPointer((realtype *)y_data,cvode->y); 1277cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 128186e87acSLisandro Dalcin 129186e87acSLisandro Dalcin for (i = 0; i < ts->max_steps; i++) { 1307cb94ee6SHong Zhang if (ts->ptime >= ts->max_time) break; 1313f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 132186e87acSLisandro Dalcin 1332bfc04deSHong Zhang if (cvode->monitorstep) { 13428adb3f7SHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP); 1352bfc04deSHong Zhang } else { 1362bfc04deSHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL); 1372bfc04deSHong Zhang } 1389f94935aSHong Zhang 1399f94935aSHong Zhang if (flag){ /* display error message */ 1409f94935aSHong Zhang switch (flag){ 1419f94935aSHong Zhang case CV_ILL_INPUT: 1429f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT"); 1439f94935aSHong Zhang break; 1449f94935aSHong Zhang case CV_TOO_CLOSE: 1459f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE"); 1469f94935aSHong Zhang break; 1473c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1489f94935aSHong Zhang PetscReal tcur; 1499f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 1509f94935aSHong Zhang ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr); 1519f94935aSHong 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); 1523c7fefeeSJed Brown } break; 1539f94935aSHong Zhang case CV_TOO_MUCH_ACC: 1549f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC"); 1559f94935aSHong Zhang break; 1569f94935aSHong Zhang case CV_ERR_FAILURE: 1579f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE"); 1589f94935aSHong Zhang break; 1599f94935aSHong Zhang case CV_CONV_FAILURE: 1609f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE"); 1619f94935aSHong Zhang break; 1629f94935aSHong Zhang case CV_LINIT_FAIL: 1639f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL"); 1649f94935aSHong Zhang break; 1659f94935aSHong Zhang case CV_LSETUP_FAIL: 1669f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL"); 1679f94935aSHong Zhang break; 1689f94935aSHong Zhang case CV_LSOLVE_FAIL: 1699f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL"); 1709f94935aSHong Zhang break; 1719f94935aSHong Zhang case CV_RHSFUNC_FAIL: 1729f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL"); 1739f94935aSHong Zhang break; 1749f94935aSHong Zhang case CV_FIRST_RHSFUNC_ERR: 1759f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR"); 1769f94935aSHong Zhang break; 1779f94935aSHong Zhang case CV_REPTD_RHSFUNC_ERR: 1789f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR"); 1799f94935aSHong Zhang break; 1809f94935aSHong Zhang case CV_UNREC_RHSFUNC_ERR: 1819f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR"); 1829f94935aSHong Zhang break; 1839f94935aSHong Zhang case CV_RTFUNC_FAIL: 1849f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL"); 1859f94935aSHong Zhang break; 1869f94935aSHong Zhang default: 1879f94935aSHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag); 1889f94935aSHong Zhang } 1899f94935aSHong Zhang } 1909f94935aSHong Zhang 1917cb94ee6SHong Zhang if (t > ts->max_time && cvode->exact_final_time) { 1927cb94ee6SHong Zhang /* interpolate to final requested time */ 1937cb94ee6SHong Zhang ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr); 1947cb94ee6SHong Zhang t = tout; 1957cb94ee6SHong Zhang } 1967cb94ee6SHong Zhang 1977cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 1987cb94ee6SHong Zhang ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr); 1997cb94ee6SHong Zhang ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr); 2007cb94ee6SHong Zhang ierr = VecResetArray(cvode->w1); CHKERRQ(ierr); 2017cb94ee6SHong Zhang ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr); 2027cb94ee6SHong Zhang ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr); 2034ac7836bSHong Zhang ierr = CVSpilsGetNumLinIters(mem, &its); 204186e87acSLisandro Dalcin ts->nonlinear_its = its; ts->linear_its = its; 205186e87acSLisandro Dalcin 206186e87acSLisandro Dalcin ts->time_step = t - ts->ptime; 207186e87acSLisandro Dalcin ts->ptime = t; 2087cb94ee6SHong Zhang ts->steps++; 209186e87acSLisandro Dalcin 2103f2090d5SJed Brown ierr = TSPostStep(ts);CHKERRQ(ierr); 2117cb94ee6SHong Zhang ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr); 2127cb94ee6SHong Zhang } 2139f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 2147cb94ee6SHong Zhang PetscFunctionReturn(0); 2157cb94ee6SHong Zhang } 2167cb94ee6SHong Zhang 2177cb94ee6SHong Zhang #undef __FUNCT__ 218277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Sundials" 219277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts) 220277b19d0SLisandro Dalcin { 221277b19d0SLisandro Dalcin TS_Sundials *cvode = (TS_Sundials*)ts->data; 222277b19d0SLisandro Dalcin PetscErrorCode ierr; 223277b19d0SLisandro Dalcin 224277b19d0SLisandro Dalcin PetscFunctionBegin; 2255c6b4a3dSJed Brown ierr = VecDestroy(&cvode->update);CHKERRQ(ierr); 226*0679a0aeSJed Brown ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr); 2275c6b4a3dSJed Brown ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr); 2285c6b4a3dSJed Brown ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr); 229277b19d0SLisandro Dalcin if (cvode->mem) {CVodeFree(&cvode->mem);} 230277b19d0SLisandro Dalcin PetscFunctionReturn(0); 231277b19d0SLisandro Dalcin } 232277b19d0SLisandro Dalcin 233277b19d0SLisandro Dalcin #undef __FUNCT__ 2347cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials" 2357cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts) 2367cb94ee6SHong Zhang { 2377cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2387cb94ee6SHong Zhang PetscErrorCode ierr; 2397cb94ee6SHong Zhang 2407cb94ee6SHong Zhang PetscFunctionBegin; 241277b19d0SLisandro Dalcin ierr = TSReset_Sundials(ts);CHKERRQ(ierr); 242a07356b8SSatish Balay ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr); 243277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 244335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","",PETSC_NULL);CHKERRQ(ierr); 245335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C","",PETSC_NULL);CHKERRQ(ierr); 246335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",PETSC_NULL);CHKERRQ(ierr); 247335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",PETSC_NULL);CHKERRQ(ierr); 248335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C","",PETSC_NULL);CHKERRQ(ierr); 249335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",PETSC_NULL);CHKERRQ(ierr); 250335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",PETSC_NULL);CHKERRQ(ierr); 251335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C","",PETSC_NULL);CHKERRQ(ierr); 252335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C","",PETSC_NULL);CHKERRQ(ierr); 253335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C","",PETSC_NULL);CHKERRQ(ierr); 254335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",PETSC_NULL);CHKERRQ(ierr); 2557cb94ee6SHong Zhang PetscFunctionReturn(0); 2567cb94ee6SHong Zhang } 2577cb94ee6SHong Zhang 2587cb94ee6SHong Zhang #undef __FUNCT__ 259214bc6a2SJed Brown #define __FUNCT__ "TSSetUp_Sundials" 260214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts) 2617cb94ee6SHong Zhang { 2627cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2637cb94ee6SHong Zhang PetscErrorCode ierr; 26416016685SHong Zhang PetscInt glosize,locsize,i,flag; 2657cb94ee6SHong Zhang PetscScalar *y_data,*parray; 26616016685SHong Zhang void *mem; 267dcbc6d53SSean Farley PC pc; 26816016685SHong Zhang const PCType pctype; 269ace3abfcSBarry Smith PetscBool pcnone; 27016016685SHong Zhang Vec sol = ts->vec_sol; 2717cb94ee6SHong Zhang 2727cb94ee6SHong Zhang PetscFunctionBegin; 2739f94935aSHong Zhang 2747cb94ee6SHong Zhang /* get the vector size */ 2757cb94ee6SHong Zhang ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr); 2767cb94ee6SHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr); 2777cb94ee6SHong Zhang 2787cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 279a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 280e32f2f54SBarry Smith if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated"); 2817cb94ee6SHong Zhang 28228adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 2837cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr); 2847cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y); 2857cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 2867cb94ee6SHong Zhang /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/ 2877cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 2887cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr); 289*0679a0aeSJed Brown ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr); 2907cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr); 291*0679a0aeSJed Brown ierr = PetscLogObjectParent(ts,cvode->ydot);CHKERRQ(ierr); 2927cb94ee6SHong Zhang 2937cb94ee6SHong Zhang /* 2947cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 2957cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 2967cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 2977cb94ee6SHong Zhang */ 2987adad957SLisandro Dalcin ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr); 2997adad957SLisandro Dalcin ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr); 3007cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr); 3017cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr); 30216016685SHong Zhang 30316016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 30416016685SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 305e32f2f54SBarry Smith if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails"); 30616016685SHong Zhang cvode->mem = mem; 30716016685SHong Zhang 30816016685SHong Zhang /* Set the pointer to user-defined data */ 30916016685SHong Zhang flag = CVodeSetUserData(mem, ts); 310e32f2f54SBarry Smith if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails"); 31116016685SHong Zhang 312fc6b9e64SJed Brown /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 313fc6b9e64SJed Brown flag = CVodeSetInitStep(mem,(realtype)ts->initial_time_step); 314fc6b9e64SJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed"); 315f1cd61daSJed Brown if (cvode->mindt > 0) { 316f1cd61daSJed Brown flag = CVodeSetMinStep(mem,(realtype)cvode->mindt); 3179f94935aSHong Zhang if (flag){ 3189f94935aSHong Zhang if (flag == CV_MEM_NULL){ 3199f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 3209f94935aSHong Zhang } else if (flag == CV_ILL_INPUT){ 3219f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size"); 3229f94935aSHong Zhang } else { 3239f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed"); 3249f94935aSHong Zhang } 3259f94935aSHong Zhang } 326f1cd61daSJed Brown } 327f1cd61daSJed Brown if (cvode->maxdt > 0) { 328f1cd61daSJed Brown flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt); 329f1cd61daSJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed"); 330f1cd61daSJed Brown } 331f1cd61daSJed Brown 33216016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 33316016685SHong Zhang * user's right hand side function in u'=f(t,u), the inital time T0, and 33416016685SHong Zhang * the initial dependent variable vector cvode->y */ 33516016685SHong Zhang flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 33616016685SHong Zhang if (flag){ 337e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag); 33816016685SHong Zhang } 33916016685SHong Zhang 3409f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 34116016685SHong Zhang flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 34216016685SHong Zhang if (flag){ 343e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag); 34416016685SHong Zhang } 34516016685SHong Zhang 3469f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 3479f94935aSHong Zhang flag = CVodeSetMaxNumSteps(mem,ts->max_steps); 34816016685SHong Zhang ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 34916016685SHong Zhang 35016016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 35116016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 352dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr); 353dcbc6d53SSean Farley ierr = PCGetType(pc,&pctype);CHKERRQ(ierr); 354dcbc6d53SSean Farley ierr = PetscTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr); 35516016685SHong Zhang if (pcnone){ 35616016685SHong Zhang flag = CVSpgmr(mem,PREC_NONE,0); 357e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 35816016685SHong Zhang } else { 35916016685SHong Zhang flag = CVSpgmr(mem,PREC_LEFT,0); 360e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 36116016685SHong Zhang 36216016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 36316016685SHong Zhang and the pointer to the user-defined block data */ 36416016685SHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 365e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 36616016685SHong Zhang } 36716016685SHong Zhang 36816016685SHong Zhang flag = CVSpilsSetGSType(mem, MODIFIED_GS); 369e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag); 3707cb94ee6SHong Zhang PetscFunctionReturn(0); 3717cb94ee6SHong Zhang } 3727cb94ee6SHong Zhang 3736fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 374dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0}; 3756fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 376dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0}; 377a04cf4d8SBarry Smith 3787cb94ee6SHong Zhang #undef __FUNCT__ 379214bc6a2SJed Brown #define __FUNCT__ "TSSetFromOptions_Sundials" 380214bc6a2SJed Brown PetscErrorCode TSSetFromOptions_Sundials(TS ts) 3817cb94ee6SHong Zhang { 3827cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 3837cb94ee6SHong Zhang PetscErrorCode ierr; 3847cb94ee6SHong Zhang int indx; 385ace3abfcSBarry Smith PetscBool flag; 3867cb94ee6SHong Zhang 3877cb94ee6SHong Zhang PetscFunctionBegin; 3887cb94ee6SHong Zhang ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr); 3896fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr); 3907cb94ee6SHong Zhang if (flag) { 3916fadb2cdSHong Zhang ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr); 3927cb94ee6SHong Zhang } 3936fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr); 3947cb94ee6SHong Zhang if (flag) { 3957cb94ee6SHong Zhang ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 3967cb94ee6SHong Zhang } 3977cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr); 3987cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr); 399f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr); 400f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr); 4017cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr); 4027cb94ee6SHong Zhang ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr); 403acfcf0e5SJed 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); 404acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr); 4057cb94ee6SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4067cb94ee6SHong Zhang PetscFunctionReturn(0); 4077cb94ee6SHong Zhang } 4087cb94ee6SHong Zhang 4097cb94ee6SHong Zhang #undef __FUNCT__ 4107cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials" 4117cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 4127cb94ee6SHong Zhang { 4137cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4147cb94ee6SHong Zhang PetscErrorCode ierr; 4157cb94ee6SHong Zhang char *type; 4167cb94ee6SHong Zhang char atype[] = "Adams"; 4177cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 418ace3abfcSBarry Smith PetscBool iascii,isstring; 4192c823083SHong Zhang long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 4202c823083SHong Zhang PetscInt qlast,qcur; 4215d47aee6SHong Zhang PetscReal hinused,hlast,hcur,tcur,tolsfac; 4227cb94ee6SHong Zhang 4237cb94ee6SHong Zhang PetscFunctionBegin; 4247cb94ee6SHong Zhang if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;} 4257cb94ee6SHong Zhang else {type = btype;} 4267cb94ee6SHong Zhang 4272692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4282692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 4297cb94ee6SHong Zhang if (iascii) { 4307cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr); 4317cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr); 4327cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr); 4337cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr); 4347cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr); 4357cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 4367cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4377cb94ee6SHong Zhang } else { 4387cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4397cb94ee6SHong Zhang } 440f1cd61daSJed Brown if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);} 441f1cd61daSJed Brown if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);} 4422c823083SHong Zhang 4435d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4445d47aee6SHong Zhang ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr); 4455d47aee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr); 4462c823083SHong Zhang ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals, 4472c823083SHong Zhang &nlinsetups,&nfails,&qlast,&qcur, 4482c823083SHong Zhang &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr); 4492c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr); 4502c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr); 4512c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr); 4522c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr); 4532c823083SHong Zhang 4542c823083SHong Zhang ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr); 4552c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr); 4562c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr); 4572c823083SHong Zhang 4582c823083SHong Zhang ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */ 4592c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr); 4602c823083SHong Zhang ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr); 4612c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr); 4622c823083SHong Zhang ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4632c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr); 4642c823083SHong Zhang ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr); 4652c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr); 4662c823083SHong Zhang ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4672c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr); 4682c823083SHong Zhang ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4692c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr); 4707cb94ee6SHong Zhang } else if (isstring) { 4717cb94ee6SHong Zhang ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 4727cb94ee6SHong Zhang } else { 473e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name); 4747cb94ee6SHong Zhang } 4757cb94ee6SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 4767cb94ee6SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 4777cb94ee6SHong Zhang PetscFunctionReturn(0); 4787cb94ee6SHong Zhang } 4797cb94ee6SHong Zhang 4807cb94ee6SHong Zhang 4817cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 4827cb94ee6SHong Zhang EXTERN_C_BEGIN 4837cb94ee6SHong Zhang #undef __FUNCT__ 4847cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials" 4857087cfbeSBarry Smith PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 4867cb94ee6SHong Zhang { 4877cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4887cb94ee6SHong Zhang 4897cb94ee6SHong Zhang PetscFunctionBegin; 4907cb94ee6SHong Zhang cvode->cvode_type = type; 4917cb94ee6SHong Zhang PetscFunctionReturn(0); 4927cb94ee6SHong Zhang } 4937cb94ee6SHong Zhang EXTERN_C_END 4947cb94ee6SHong Zhang 4957cb94ee6SHong Zhang EXTERN_C_BEGIN 4967cb94ee6SHong Zhang #undef __FUNCT__ 4977cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials" 4987087cfbeSBarry Smith PetscErrorCode TSSundialsSetGMRESRestart_Sundials(TS ts,int restart) 4997cb94ee6SHong Zhang { 5007cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5017cb94ee6SHong Zhang 5027cb94ee6SHong Zhang PetscFunctionBegin; 5037cb94ee6SHong Zhang cvode->restart = restart; 5047cb94ee6SHong Zhang PetscFunctionReturn(0); 5057cb94ee6SHong Zhang } 5067cb94ee6SHong Zhang EXTERN_C_END 5077cb94ee6SHong Zhang 5087cb94ee6SHong Zhang EXTERN_C_BEGIN 5097cb94ee6SHong Zhang #undef __FUNCT__ 5107cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials" 5117087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,double tol) 5127cb94ee6SHong Zhang { 5137cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5147cb94ee6SHong Zhang 5157cb94ee6SHong Zhang PetscFunctionBegin; 5167cb94ee6SHong Zhang cvode->linear_tol = tol; 5177cb94ee6SHong Zhang PetscFunctionReturn(0); 5187cb94ee6SHong Zhang } 5197cb94ee6SHong Zhang EXTERN_C_END 5207cb94ee6SHong Zhang 5217cb94ee6SHong Zhang EXTERN_C_BEGIN 5227cb94ee6SHong Zhang #undef __FUNCT__ 5237cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials" 5247087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 5257cb94ee6SHong Zhang { 5267cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5277cb94ee6SHong Zhang 5287cb94ee6SHong Zhang PetscFunctionBegin; 5297cb94ee6SHong Zhang cvode->gtype = type; 5307cb94ee6SHong Zhang PetscFunctionReturn(0); 5317cb94ee6SHong Zhang } 5327cb94ee6SHong Zhang EXTERN_C_END 5337cb94ee6SHong Zhang 5347cb94ee6SHong Zhang EXTERN_C_BEGIN 5357cb94ee6SHong Zhang #undef __FUNCT__ 5367cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials" 5377087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel) 5387cb94ee6SHong Zhang { 5397cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5407cb94ee6SHong Zhang 5417cb94ee6SHong Zhang PetscFunctionBegin; 5427cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 5437cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 5447cb94ee6SHong Zhang PetscFunctionReturn(0); 5457cb94ee6SHong Zhang } 5467cb94ee6SHong Zhang EXTERN_C_END 5477cb94ee6SHong Zhang 5487cb94ee6SHong Zhang EXTERN_C_BEGIN 5497cb94ee6SHong Zhang #undef __FUNCT__ 550f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials" 5517087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 552f1cd61daSJed Brown { 553f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 554f1cd61daSJed Brown 555f1cd61daSJed Brown PetscFunctionBegin; 556f1cd61daSJed Brown cvode->mindt = mindt; 557f1cd61daSJed Brown PetscFunctionReturn(0); 558f1cd61daSJed Brown } 559f1cd61daSJed Brown EXTERN_C_END 560f1cd61daSJed Brown 561f1cd61daSJed Brown EXTERN_C_BEGIN 562f1cd61daSJed Brown #undef __FUNCT__ 563f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials" 5647087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 565f1cd61daSJed Brown { 566f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 567f1cd61daSJed Brown 568f1cd61daSJed Brown PetscFunctionBegin; 569f1cd61daSJed Brown cvode->maxdt = maxdt; 570f1cd61daSJed Brown PetscFunctionReturn(0); 571f1cd61daSJed Brown } 572f1cd61daSJed Brown EXTERN_C_END 573f1cd61daSJed Brown 574f1cd61daSJed Brown EXTERN_C_BEGIN 575f1cd61daSJed Brown #undef __FUNCT__ 5767cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials" 5777087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc) 5787cb94ee6SHong Zhang { 579dcbc6d53SSean Farley SNES snes; 580dcbc6d53SSean Farley KSP ksp; 581dcbc6d53SSean Farley PetscErrorCode ierr; 5827cb94ee6SHong Zhang 5837cb94ee6SHong Zhang PetscFunctionBegin; 584dcbc6d53SSean Farley ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 585dcbc6d53SSean Farley ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 586dcbc6d53SSean Farley ierr = KSPGetPC(ksp,pc); CHKERRQ(ierr); 5877cb94ee6SHong Zhang PetscFunctionReturn(0); 5887cb94ee6SHong Zhang } 5897cb94ee6SHong Zhang EXTERN_C_END 5907cb94ee6SHong Zhang 5917cb94ee6SHong Zhang EXTERN_C_BEGIN 5927cb94ee6SHong Zhang #undef __FUNCT__ 5937cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials" 5947087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 5957cb94ee6SHong Zhang { 5967cb94ee6SHong Zhang PetscFunctionBegin; 5972c823083SHong Zhang if (nonlin) *nonlin = ts->nonlinear_its; 5982c823083SHong Zhang if (lin) *lin = ts->linear_its; 5997cb94ee6SHong Zhang PetscFunctionReturn(0); 6007cb94ee6SHong Zhang } 6017cb94ee6SHong Zhang EXTERN_C_END 6027cb94ee6SHong Zhang 6037cb94ee6SHong Zhang EXTERN_C_BEGIN 6047cb94ee6SHong Zhang #undef __FUNCT__ 6057cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials" 6067087cfbeSBarry Smith PetscErrorCode TSSundialsSetExactFinalTime_Sundials(TS ts,PetscBool s) 6077cb94ee6SHong Zhang { 6087cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 6097cb94ee6SHong Zhang 6107cb94ee6SHong Zhang PetscFunctionBegin; 6117cb94ee6SHong Zhang cvode->exact_final_time = s; 6127cb94ee6SHong Zhang PetscFunctionReturn(0); 6137cb94ee6SHong Zhang } 6147cb94ee6SHong Zhang EXTERN_C_END 6152bfc04deSHong Zhang 6162bfc04deSHong Zhang EXTERN_C_BEGIN 6172bfc04deSHong Zhang #undef __FUNCT__ 6182bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials" 6197087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s) 6202bfc04deSHong Zhang { 6212bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 6222bfc04deSHong Zhang 6232bfc04deSHong Zhang PetscFunctionBegin; 6242bfc04deSHong Zhang cvode->monitorstep = s; 6252bfc04deSHong Zhang PetscFunctionReturn(0); 6262bfc04deSHong Zhang } 6272bfc04deSHong Zhang EXTERN_C_END 6287cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 6297cb94ee6SHong Zhang 6307cb94ee6SHong Zhang #undef __FUNCT__ 6317cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations" 6327cb94ee6SHong Zhang /*@C 6337cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 6347cb94ee6SHong Zhang 6357cb94ee6SHong Zhang Not Collective 6367cb94ee6SHong Zhang 6377cb94ee6SHong Zhang Input parameters: 6387cb94ee6SHong Zhang . ts - the time-step context 6397cb94ee6SHong Zhang 6407cb94ee6SHong Zhang Output Parameters: 6417cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 6427cb94ee6SHong Zhang - lin - number of linear iterations 6437cb94ee6SHong Zhang 6447cb94ee6SHong Zhang Level: advanced 6457cb94ee6SHong Zhang 6467cb94ee6SHong Zhang Notes: 6477cb94ee6SHong Zhang These return the number since the creation of the TS object 6487cb94ee6SHong Zhang 6497cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations 6507cb94ee6SHong Zhang 6517cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6527cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 6537cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 654f1cd61daSJed Brown TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSundialsSetExactFinalTime() 6557cb94ee6SHong Zhang 6567cb94ee6SHong Zhang @*/ 6577087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 6587cb94ee6SHong Zhang { 6594ac538c5SBarry Smith PetscErrorCode ierr; 6607cb94ee6SHong Zhang 6617cb94ee6SHong Zhang PetscFunctionBegin; 6624ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr); 6637cb94ee6SHong Zhang PetscFunctionReturn(0); 6647cb94ee6SHong Zhang } 6657cb94ee6SHong Zhang 6667cb94ee6SHong Zhang #undef __FUNCT__ 6677cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType" 6687cb94ee6SHong Zhang /*@ 6697cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 6707cb94ee6SHong Zhang 6713f9fe445SBarry Smith Logically Collective on TS 6727cb94ee6SHong Zhang 6737cb94ee6SHong Zhang Input parameters: 6747cb94ee6SHong Zhang + ts - the time-step context 6757cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 6767cb94ee6SHong Zhang 6777cb94ee6SHong Zhang Level: intermediate 6787cb94ee6SHong Zhang 6797cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula 6807cb94ee6SHong Zhang 6817cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetGMRESRestart(), 6827cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 6837cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 6847cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 6857cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 6867cb94ee6SHong Zhang @*/ 6877087cfbeSBarry Smith PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type) 6887cb94ee6SHong Zhang { 6894ac538c5SBarry Smith PetscErrorCode ierr; 6907cb94ee6SHong Zhang 6917cb94ee6SHong Zhang PetscFunctionBegin; 6924ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr); 6937cb94ee6SHong Zhang PetscFunctionReturn(0); 6947cb94ee6SHong Zhang } 6957cb94ee6SHong Zhang 6967cb94ee6SHong Zhang #undef __FUNCT__ 6977cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart" 6987cb94ee6SHong Zhang /*@ 6997cb94ee6SHong Zhang TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by 7007cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 7017cb94ee6SHong Zhang this is ALSO the maximum number of GMRES steps that will be used. 7027cb94ee6SHong Zhang 7033f9fe445SBarry Smith Logically Collective on TS 7047cb94ee6SHong Zhang 7057cb94ee6SHong Zhang Input parameters: 7067cb94ee6SHong Zhang + ts - the time-step context 7077cb94ee6SHong Zhang - restart - number of direction vectors (the restart size). 7087cb94ee6SHong Zhang 7097cb94ee6SHong Zhang Level: advanced 7107cb94ee6SHong Zhang 7117cb94ee6SHong Zhang .keywords: GMRES, restart 7127cb94ee6SHong Zhang 7137cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 7147cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 7157cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7167cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 7177cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 7187cb94ee6SHong Zhang 7197cb94ee6SHong Zhang @*/ 7207087cfbeSBarry Smith PetscErrorCode TSSundialsSetGMRESRestart(TS ts,int restart) 7217cb94ee6SHong Zhang { 7224ac538c5SBarry Smith PetscErrorCode ierr; 7237cb94ee6SHong Zhang 7247cb94ee6SHong Zhang PetscFunctionBegin; 725c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(ts,restart,2); 7264ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetGMRESRestart_C",(TS,int),(ts,restart));CHKERRQ(ierr); 7277cb94ee6SHong Zhang PetscFunctionReturn(0); 7287cb94ee6SHong Zhang } 7297cb94ee6SHong Zhang 7307cb94ee6SHong Zhang #undef __FUNCT__ 7317cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance" 7327cb94ee6SHong Zhang /*@ 7337cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 7347cb94ee6SHong Zhang system by SUNDIALS. 7357cb94ee6SHong Zhang 7363f9fe445SBarry Smith Logically Collective on TS 7377cb94ee6SHong Zhang 7387cb94ee6SHong Zhang Input parameters: 7397cb94ee6SHong Zhang + ts - the time-step context 7407cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 7417cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 7427cb94ee6SHong Zhang 7437cb94ee6SHong Zhang Level: advanced 7447cb94ee6SHong Zhang 7457cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS 7467cb94ee6SHong Zhang 7477cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7487cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 7497cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7507cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 7517cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 7527cb94ee6SHong Zhang 7537cb94ee6SHong Zhang @*/ 7547087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance(TS ts,double tol) 7557cb94ee6SHong Zhang { 7564ac538c5SBarry Smith PetscErrorCode ierr; 7577cb94ee6SHong Zhang 7587cb94ee6SHong Zhang PetscFunctionBegin; 759c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,tol,2); 7604ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr); 7617cb94ee6SHong Zhang PetscFunctionReturn(0); 7627cb94ee6SHong Zhang } 7637cb94ee6SHong Zhang 7647cb94ee6SHong Zhang #undef __FUNCT__ 7657cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType" 7667cb94ee6SHong Zhang /*@ 7677cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 7687cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 7697cb94ee6SHong Zhang 7703f9fe445SBarry Smith Logically Collective on TS 7717cb94ee6SHong Zhang 7727cb94ee6SHong Zhang Input parameters: 7737cb94ee6SHong Zhang + ts - the time-step context 7747cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 7757cb94ee6SHong Zhang 7767cb94ee6SHong Zhang Level: advanced 7777cb94ee6SHong Zhang 7787cb94ee6SHong Zhang .keywords: Sundials, orthogonalization 7797cb94ee6SHong Zhang 7807cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7817cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 7827cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 7837cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 7847cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 7857cb94ee6SHong Zhang 7867cb94ee6SHong Zhang @*/ 7877087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 7887cb94ee6SHong Zhang { 7894ac538c5SBarry Smith PetscErrorCode ierr; 7907cb94ee6SHong Zhang 7917cb94ee6SHong Zhang PetscFunctionBegin; 7924ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr); 7937cb94ee6SHong Zhang PetscFunctionReturn(0); 7947cb94ee6SHong Zhang } 7957cb94ee6SHong Zhang 7967cb94ee6SHong Zhang #undef __FUNCT__ 7977cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance" 7987cb94ee6SHong Zhang /*@ 7997cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 8007cb94ee6SHong Zhang Sundials for error control. 8017cb94ee6SHong Zhang 8023f9fe445SBarry Smith Logically Collective on TS 8037cb94ee6SHong Zhang 8047cb94ee6SHong Zhang Input parameters: 8057cb94ee6SHong Zhang + ts - the time-step context 8067cb94ee6SHong Zhang . aabs - the absolute tolerance 8077cb94ee6SHong Zhang - rel - the relative tolerance 8087cb94ee6SHong Zhang 8097cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 8107cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 8117cb94ee6SHong Zhang 8127cb94ee6SHong Zhang Level: intermediate 8137cb94ee6SHong Zhang 8147cb94ee6SHong Zhang .keywords: Sundials, tolerance 8157cb94ee6SHong Zhang 8167cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8177cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 8187cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8197cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 8207cb94ee6SHong Zhang TSSundialsSetExactFinalTime() 8217cb94ee6SHong Zhang 8227cb94ee6SHong Zhang @*/ 8237087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance(TS ts,double aabs,double rel) 8247cb94ee6SHong Zhang { 8254ac538c5SBarry Smith PetscErrorCode ierr; 8267cb94ee6SHong Zhang 8277cb94ee6SHong Zhang PetscFunctionBegin; 8284ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr); 8297cb94ee6SHong Zhang PetscFunctionReturn(0); 8307cb94ee6SHong Zhang } 8317cb94ee6SHong Zhang 8327cb94ee6SHong Zhang #undef __FUNCT__ 8337cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC" 8347cb94ee6SHong Zhang /*@ 8357cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 8367cb94ee6SHong Zhang 8377cb94ee6SHong Zhang Input Parameter: 8387cb94ee6SHong Zhang . ts - the time-step context 8397cb94ee6SHong Zhang 8407cb94ee6SHong Zhang Output Parameter: 8417cb94ee6SHong Zhang . pc - the preconditioner context 8427cb94ee6SHong Zhang 8437cb94ee6SHong Zhang Level: advanced 8447cb94ee6SHong Zhang 8457cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8467cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 8477cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8487cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 8497cb94ee6SHong Zhang @*/ 8507087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC(TS ts,PC *pc) 8517cb94ee6SHong Zhang { 8524ac538c5SBarry Smith PetscErrorCode ierr; 8537cb94ee6SHong Zhang 8547cb94ee6SHong Zhang PetscFunctionBegin; 8554ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr); 8567cb94ee6SHong Zhang PetscFunctionReturn(0); 8577cb94ee6SHong Zhang } 8587cb94ee6SHong Zhang 8597cb94ee6SHong Zhang #undef __FUNCT__ 8607cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime" 8617cb94ee6SHong Zhang /*@ 8627cb94ee6SHong Zhang TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the 8637cb94ee6SHong Zhang exact final time requested by the user or just returns it at the final time 8647cb94ee6SHong Zhang it computed. (Defaults to true). 8657cb94ee6SHong Zhang 8667cb94ee6SHong Zhang Input Parameter: 8677cb94ee6SHong Zhang + ts - the time-step context 8687cb94ee6SHong Zhang - ft - PETSC_TRUE if interpolates, else PETSC_FALSE 8697cb94ee6SHong Zhang 8707cb94ee6SHong Zhang Level: beginner 8717cb94ee6SHong Zhang 8727cb94ee6SHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8737cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 8747cb94ee6SHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 8757cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 8767cb94ee6SHong Zhang @*/ 8777087cfbeSBarry Smith PetscErrorCode TSSundialsSetExactFinalTime(TS ts,PetscBool ft) 8787cb94ee6SHong Zhang { 8794ac538c5SBarry Smith PetscErrorCode ierr; 8807cb94ee6SHong Zhang 8817cb94ee6SHong Zhang PetscFunctionBegin; 8824ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetExactFinalTime_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 8837cb94ee6SHong Zhang PetscFunctionReturn(0); 8847cb94ee6SHong Zhang } 8857cb94ee6SHong Zhang 8862bfc04deSHong Zhang #undef __FUNCT__ 887f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep" 888f1cd61daSJed Brown /*@ 889f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 890f1cd61daSJed Brown 891f1cd61daSJed Brown Input Parameter: 892f1cd61daSJed Brown + ts - the time-step context 893f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 894f1cd61daSJed Brown 895fc6b9e64SJed Brown Note: 896fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 897fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 898fc6b9e64SJed Brown 899f1cd61daSJed Brown Level: beginner 900f1cd61daSJed Brown 901f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 902f1cd61daSJed Brown @*/ 9037087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 904f1cd61daSJed Brown { 9054ac538c5SBarry Smith PetscErrorCode ierr; 906f1cd61daSJed Brown 907f1cd61daSJed Brown PetscFunctionBegin; 9084ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr); 909f1cd61daSJed Brown PetscFunctionReturn(0); 910f1cd61daSJed Brown } 911f1cd61daSJed Brown 912f1cd61daSJed Brown #undef __FUNCT__ 913f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep" 914f1cd61daSJed Brown /*@ 915f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 916f1cd61daSJed Brown 917f1cd61daSJed Brown Input Parameter: 918f1cd61daSJed Brown + ts - the time-step context 919f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 920f1cd61daSJed Brown 921f1cd61daSJed Brown Level: beginner 922f1cd61daSJed Brown 923f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 924f1cd61daSJed Brown @*/ 9257087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 926f1cd61daSJed Brown { 9274ac538c5SBarry Smith PetscErrorCode ierr; 928f1cd61daSJed Brown 929f1cd61daSJed Brown PetscFunctionBegin; 9304ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 931f1cd61daSJed Brown PetscFunctionReturn(0); 932f1cd61daSJed Brown } 933f1cd61daSJed Brown 934f1cd61daSJed Brown #undef __FUNCT__ 9352bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps" 9362bfc04deSHong Zhang /*@ 9372bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 9382bfc04deSHong Zhang 9392bfc04deSHong Zhang Input Parameter: 9402bfc04deSHong Zhang + ts - the time-step context 9412bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 9422bfc04deSHong Zhang 9432bfc04deSHong Zhang Level: beginner 9442bfc04deSHong Zhang 9452bfc04deSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 9462bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 9472bfc04deSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), 9482bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 9492bfc04deSHong Zhang @*/ 9507087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft) 9512bfc04deSHong Zhang { 9524ac538c5SBarry Smith PetscErrorCode ierr; 9532bfc04deSHong Zhang 9542bfc04deSHong Zhang PetscFunctionBegin; 9554ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 9562bfc04deSHong Zhang PetscFunctionReturn(0); 9572bfc04deSHong Zhang } 9587cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 9597cb94ee6SHong Zhang /*MC 96096f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 9617cb94ee6SHong Zhang 9627cb94ee6SHong Zhang Options Database: 9637cb94ee6SHong Zhang + -ts_sundials_type <bdf,adams> 9647cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 9657cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 9667cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 9677cb94ee6SHong Zhang . -ts_sundials_linear_tolerance <tol> 9687cb94ee6SHong Zhang . -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions 96916016685SHong Zhang . -ts_sundials_exact_final_time - Interpolate output to stop exactly at the final time 97016016685SHong Zhang - -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps 97116016685SHong Zhang 9727cb94ee6SHong Zhang 9737cb94ee6SHong Zhang Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply 9747cb94ee6SHong Zhang only PETSc PC options 9757cb94ee6SHong Zhang 9767cb94ee6SHong Zhang Level: beginner 9777cb94ee6SHong Zhang 9787cb94ee6SHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(), 9797cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime() 9807cb94ee6SHong Zhang 9817cb94ee6SHong Zhang M*/ 9827cb94ee6SHong Zhang EXTERN_C_BEGIN 9837cb94ee6SHong Zhang #undef __FUNCT__ 9847cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials" 9857087cfbeSBarry Smith PetscErrorCode TSCreate_Sundials(TS ts) 9867cb94ee6SHong Zhang { 9877cb94ee6SHong Zhang TS_Sundials *cvode; 9887cb94ee6SHong Zhang PetscErrorCode ierr; 9897cb94ee6SHong Zhang 9907cb94ee6SHong Zhang PetscFunctionBegin; 991277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Sundials; 99228adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 99328adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 994214bc6a2SJed Brown ts->ops->setup = TSSetUp_Sundials; 995214bc6a2SJed Brown ts->ops->step = TSStep_Sundials; 996214bc6a2SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Sundials; 9977cb94ee6SHong Zhang 99838f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr); 9997cb94ee6SHong Zhang ts->data = (void*)cvode; 10006fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 10016fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 10027cb94ee6SHong Zhang cvode->restart = 5; 10037cb94ee6SHong Zhang cvode->linear_tol = .05; 10047cb94ee6SHong Zhang 10052bfc04deSHong Zhang cvode->exact_final_time = PETSC_TRUE; 10062bfc04deSHong Zhang cvode->monitorstep = PETSC_FALSE; 10077cb94ee6SHong Zhang 1008a07356b8SSatish Balay ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr); 1009f1cd61daSJed Brown 1010f1cd61daSJed Brown cvode->mindt = -1.; 1011f1cd61daSJed Brown cvode->maxdt = -1.; 1012f1cd61daSJed Brown 10137cb94ee6SHong Zhang /* set tolerance for Sundials */ 10147cb94ee6SHong Zhang cvode->reltol = 1e-6; 10152c823083SHong Zhang cvode->abstol = 1e-6; 10167cb94ee6SHong Zhang 10177cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials", 10187cb94ee6SHong Zhang TSSundialsSetType_Sundials);CHKERRQ(ierr); 10197cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C", 10207cb94ee6SHong Zhang "TSSundialsSetGMRESRestart_Sundials", 10217cb94ee6SHong Zhang TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr); 10227cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C", 10237cb94ee6SHong Zhang "TSSundialsSetLinearTolerance_Sundials", 10247cb94ee6SHong Zhang TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 10257cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C", 10267cb94ee6SHong Zhang "TSSundialsSetGramSchmidtType_Sundials", 10277cb94ee6SHong Zhang TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 10287cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C", 10297cb94ee6SHong Zhang "TSSundialsSetTolerance_Sundials", 10307cb94ee6SHong Zhang TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 1031f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C", 1032f1cd61daSJed Brown "TSSundialsSetMinTimeStep_Sundials", 1033f1cd61daSJed Brown TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr); 1034f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C", 1035f1cd61daSJed Brown "TSSundialsSetMaxTimeStep_Sundials", 1036f1cd61daSJed Brown TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr); 10377cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C", 10387cb94ee6SHong Zhang "TSSundialsGetPC_Sundials", 10397cb94ee6SHong Zhang TSSundialsGetPC_Sundials);CHKERRQ(ierr); 10407cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C", 10417cb94ee6SHong Zhang "TSSundialsGetIterations_Sundials", 10427cb94ee6SHong Zhang TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 10437cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C", 10447cb94ee6SHong Zhang "TSSundialsSetExactFinalTime_Sundials", 10457cb94ee6SHong Zhang TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr); 10462bfc04deSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C", 10472bfc04deSHong Zhang "TSSundialsMonitorInternalSteps_Sundials", 10482bfc04deSHong Zhang TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr); 10492bfc04deSHong Zhang 10507cb94ee6SHong Zhang PetscFunctionReturn(0); 10517cb94ee6SHong Zhang } 10527cb94ee6SHong Zhang EXTERN_C_END 1053