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; 250679a0aeSJed Brown Mat J,P; 260679a0aeSJed Brown Vec yy = cvode->w1,yydot = cvode->ydot; 270679a0aeSJed Brown PetscReal gm = (PetscReal)_gamma; 287cb94ee6SHong Zhang MatStructure str = DIFFERENT_NONZERO_PATTERN; 297cb94ee6SHong Zhang PetscScalar *y_data; 307cb94ee6SHong Zhang 317cb94ee6SHong Zhang PetscFunctionBegin; 320679a0aeSJed 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); 350679a0aeSJed Brown ierr = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */ 360679a0aeSJed Brown /* compute the shifted Jacobian (1/gm)*I + Jrest */ 370679a0aeSJed Brown ierr = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,&J,&P,&str,PETSC_FALSE);CHKERRQ(ierr); 387cb94ee6SHong Zhang ierr = VecResetArray(yy); CHKERRQ(ierr); 390679a0aeSJed 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); 420679a0aeSJed 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; 845fcef5e4SJed Brown DM dm; 85*942e3340SBarry Smith DMTS tsdm; 865fcef5e4SJed Brown TSIFunction ifunction; 877adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)ts)->comm; 887cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 890679a0aeSJed Brown Vec yy = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot; 907cb94ee6SHong Zhang PetscScalar *y_data,*ydot_data; 917cb94ee6SHong Zhang PetscErrorCode ierr; 927cb94ee6SHong Zhang 937cb94ee6SHong Zhang PetscFunctionBegin; 945ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 957cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(y); 967cb94ee6SHong Zhang ydot_data = (PetscScalar *) N_VGetArrayPointer(ydot); 974ac7836bSHong Zhang ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr); 984ac7836bSHong Zhang ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr); 994ac7836bSHong Zhang 1005fcef5e4SJed Brown /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */ 1015fcef5e4SJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 102*942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 1035fcef5e4SJed Brown ierr = DMTSGetIFunction(dm,&ifunction,PETSC_NULL);CHKERRQ(ierr); 1045fcef5e4SJed Brown if (!ifunction) { 105bc0cc02bSJed Brown ierr = TSComputeRHSFunction(ts,t,yy,yyd);CHKERRQ(ierr); 106bc0cc02bSJed Brown } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */ 107bc0cc02bSJed Brown ierr = VecZeroEntries(yydot);CHKERRQ(ierr); 1080679a0aeSJed Brown ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE); CHKERRABORT(comm,ierr); 1090679a0aeSJed Brown ierr = VecScale(yyd,-1.);CHKERRQ(ierr); 110bc0cc02bSJed Brown } 1117cb94ee6SHong Zhang ierr = VecResetArray(yy); CHKERRABORT(comm,ierr); 1127cb94ee6SHong Zhang ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr); 1134ac7836bSHong Zhang PetscFunctionReturn(0); 1147cb94ee6SHong Zhang } 1157cb94ee6SHong Zhang 1167cb94ee6SHong Zhang /* 117b4eba00bSSean Farley TSStep_Sundials - Calls Sundials to integrate the ODE. 1187cb94ee6SHong Zhang */ 1197cb94ee6SHong Zhang #undef __FUNCT__ 120b4eba00bSSean Farley #define __FUNCT__ "TSStep_Sundials" 121b4eba00bSSean Farley PetscErrorCode TSStep_Sundials(TS ts) 1227cb94ee6SHong Zhang { 1237cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1247cb94ee6SHong Zhang PetscErrorCode ierr; 125b4eba00bSSean Farley PetscInt flag; 1269f94935aSHong Zhang long int its,nsteps; 1277cb94ee6SHong Zhang realtype t,tout; 1287cb94ee6SHong Zhang PetscScalar *y_data; 1297cb94ee6SHong Zhang void *mem; 1307cb94ee6SHong Zhang 1317cb94ee6SHong Zhang PetscFunctionBegin; 13216016685SHong Zhang mem = cvode->mem; 1337cb94ee6SHong Zhang tout = ts->max_time; 1347cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr); 1357cb94ee6SHong Zhang N_VSetArrayPointer((realtype *)y_data,cvode->y); 1367cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 137186e87acSLisandro Dalcin 1383f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 139186e87acSLisandro Dalcin 1402bfc04deSHong Zhang if (cvode->monitorstep) { 141b8123daeSJed Brown /* We would like to call TSPreStep() when starting each step (including rejections) and TSPreStage() before each 142b8123daeSJed Brown * stage solve, but CVode does not appear to support this. */ 14328adb3f7SHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP); 1442bfc04deSHong Zhang } else { 1452bfc04deSHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL); 1462bfc04deSHong Zhang } 1479f94935aSHong Zhang 1489f94935aSHong Zhang if (flag){ /* display error message */ 1499f94935aSHong Zhang switch (flag){ 1509f94935aSHong Zhang case CV_ILL_INPUT: 1519f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT"); 1529f94935aSHong Zhang break; 1539f94935aSHong Zhang case CV_TOO_CLOSE: 1549f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE"); 1559f94935aSHong Zhang break; 1563c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1579f94935aSHong Zhang PetscReal tcur; 1589f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 1599f94935aSHong Zhang ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr); 1609f94935aSHong 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); 1613c7fefeeSJed Brown } break; 1629f94935aSHong Zhang case CV_TOO_MUCH_ACC: 1639f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC"); 1649f94935aSHong Zhang break; 1659f94935aSHong Zhang case CV_ERR_FAILURE: 1669f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE"); 1679f94935aSHong Zhang break; 1689f94935aSHong Zhang case CV_CONV_FAILURE: 1699f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE"); 1709f94935aSHong Zhang break; 1719f94935aSHong Zhang case CV_LINIT_FAIL: 1729f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL"); 1739f94935aSHong Zhang break; 1749f94935aSHong Zhang case CV_LSETUP_FAIL: 1759f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL"); 1769f94935aSHong Zhang break; 1779f94935aSHong Zhang case CV_LSOLVE_FAIL: 1789f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL"); 1799f94935aSHong Zhang break; 1809f94935aSHong Zhang case CV_RHSFUNC_FAIL: 1819f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL"); 1829f94935aSHong Zhang break; 1839f94935aSHong Zhang case CV_FIRST_RHSFUNC_ERR: 1849f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR"); 1859f94935aSHong Zhang break; 1869f94935aSHong Zhang case CV_REPTD_RHSFUNC_ERR: 1879f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR"); 1889f94935aSHong Zhang break; 1899f94935aSHong Zhang case CV_UNREC_RHSFUNC_ERR: 1909f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR"); 1919f94935aSHong Zhang break; 1929f94935aSHong Zhang case CV_RTFUNC_FAIL: 1939f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL"); 1949f94935aSHong Zhang break; 1959f94935aSHong Zhang default: 1969f94935aSHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag); 1979f94935aSHong Zhang } 1989f94935aSHong Zhang } 1999f94935aSHong Zhang 2007cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 2017cb94ee6SHong Zhang ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr); 2027cb94ee6SHong Zhang ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr); 2037cb94ee6SHong Zhang ierr = VecResetArray(cvode->w1); CHKERRQ(ierr); 204bc0cc02bSJed Brown ierr = VecCopy(cvode->update,ts->vec_sol);CHKERRQ(ierr); 2057cb94ee6SHong Zhang ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr); 2064ac7836bSHong Zhang ierr = CVSpilsGetNumLinIters(mem, &its); 2075ef26d82SJed Brown ts->snes_its = its; ts->ksp_its = its; 208186e87acSLisandro Dalcin 209186e87acSLisandro Dalcin ts->time_step = t - ts->ptime; 210186e87acSLisandro Dalcin ts->ptime = t; 2117cb94ee6SHong Zhang ts->steps++; 212186e87acSLisandro Dalcin 2139f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 214b4eba00bSSean Farley if (!cvode->monitorstep) ts->steps = nsteps; 215b4eba00bSSean Farley PetscFunctionReturn(0); 216b4eba00bSSean Farley } 217b4eba00bSSean Farley 218b4eba00bSSean Farley #undef __FUNCT__ 219b4eba00bSSean Farley #define __FUNCT__ "TSInterpolate_Sundials" 220b4eba00bSSean Farley static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X) 221b4eba00bSSean Farley { 222b4eba00bSSean Farley TS_Sundials *cvode = (TS_Sundials*)ts->data; 223b4eba00bSSean Farley N_Vector y; 224b4eba00bSSean Farley PetscErrorCode ierr; 225b4eba00bSSean Farley PetscScalar *x_data; 226b4eba00bSSean Farley PetscInt glosize,locsize; 227b4eba00bSSean Farley 228b4eba00bSSean Farley PetscFunctionBegin; 229b4eba00bSSean Farley 230b4eba00bSSean Farley /* get the vector size */ 231b4eba00bSSean Farley ierr = VecGetSize(X,&glosize);CHKERRQ(ierr); 232b4eba00bSSean Farley ierr = VecGetLocalSize(X,&locsize);CHKERRQ(ierr); 233b4eba00bSSean Farley 234b4eba00bSSean Farley /* allocate the memory for N_Vec y */ 235b4eba00bSSean Farley y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 236b4eba00bSSean Farley if (!y) SETERRQ(PETSC_COMM_SELF,1,"Interpolated y is not allocated"); 237b4eba00bSSean Farley 238b4eba00bSSean Farley ierr = VecGetArray(X,&x_data);CHKERRQ(ierr); 239b4eba00bSSean Farley N_VSetArrayPointer((realtype *)x_data,y); 240b4eba00bSSean Farley ierr = CVodeGetDky(cvode->mem,t,0,y);CHKERRQ(ierr); 241b4eba00bSSean Farley ierr = VecRestoreArray(X,&x_data);CHKERRQ(ierr); 242b4eba00bSSean Farley 2437cb94ee6SHong Zhang PetscFunctionReturn(0); 2447cb94ee6SHong Zhang } 2457cb94ee6SHong Zhang 2467cb94ee6SHong Zhang #undef __FUNCT__ 247277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Sundials" 248277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts) 249277b19d0SLisandro Dalcin { 250277b19d0SLisandro Dalcin TS_Sundials *cvode = (TS_Sundials*)ts->data; 251277b19d0SLisandro Dalcin PetscErrorCode ierr; 252277b19d0SLisandro Dalcin 253277b19d0SLisandro Dalcin PetscFunctionBegin; 2545c6b4a3dSJed Brown ierr = VecDestroy(&cvode->update);CHKERRQ(ierr); 2550679a0aeSJed Brown ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr); 2565c6b4a3dSJed Brown ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr); 2575c6b4a3dSJed Brown ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr); 258277b19d0SLisandro Dalcin if (cvode->mem) {CVodeFree(&cvode->mem);} 259277b19d0SLisandro Dalcin PetscFunctionReturn(0); 260277b19d0SLisandro Dalcin } 261277b19d0SLisandro Dalcin 262277b19d0SLisandro Dalcin #undef __FUNCT__ 2637cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials" 2647cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts) 2657cb94ee6SHong Zhang { 2667cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2677cb94ee6SHong Zhang PetscErrorCode ierr; 2687cb94ee6SHong Zhang 2697cb94ee6SHong Zhang PetscFunctionBegin; 270277b19d0SLisandro Dalcin ierr = TSReset_Sundials(ts);CHKERRQ(ierr); 271a07356b8SSatish Balay ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr); 272277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 273335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","",PETSC_NULL);CHKERRQ(ierr); 274f61b2b6aSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C","",PETSC_NULL);CHKERRQ(ierr); 275335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",PETSC_NULL);CHKERRQ(ierr); 276335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",PETSC_NULL);CHKERRQ(ierr); 277335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C","",PETSC_NULL);CHKERRQ(ierr); 278335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",PETSC_NULL);CHKERRQ(ierr); 279335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",PETSC_NULL);CHKERRQ(ierr); 280335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C","",PETSC_NULL);CHKERRQ(ierr); 281335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C","",PETSC_NULL);CHKERRQ(ierr); 282335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",PETSC_NULL);CHKERRQ(ierr); 2837cb94ee6SHong Zhang PetscFunctionReturn(0); 2847cb94ee6SHong Zhang } 2857cb94ee6SHong Zhang 2867cb94ee6SHong Zhang #undef __FUNCT__ 287214bc6a2SJed Brown #define __FUNCT__ "TSSetUp_Sundials" 288214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts) 2897cb94ee6SHong Zhang { 2907cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2917cb94ee6SHong Zhang PetscErrorCode ierr; 29216016685SHong Zhang PetscInt glosize,locsize,i,flag; 2937cb94ee6SHong Zhang PetscScalar *y_data,*parray; 29416016685SHong Zhang void *mem; 295dcbc6d53SSean Farley PC pc; 29619fd82e9SBarry Smith PCType pctype; 297ace3abfcSBarry Smith PetscBool pcnone; 2987cb94ee6SHong Zhang 2997cb94ee6SHong Zhang PetscFunctionBegin; 3007cb94ee6SHong Zhang /* get the vector size */ 3017cb94ee6SHong Zhang ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr); 3027cb94ee6SHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr); 3037cb94ee6SHong Zhang 3047cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 305a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 306e32f2f54SBarry Smith if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated"); 3077cb94ee6SHong Zhang 30828adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 3097cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr); 3107cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y); 3117cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 3127cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 31342528757SHong Zhang 3147cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr); 3150679a0aeSJed Brown ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr); 3167cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr); 3170679a0aeSJed Brown ierr = PetscLogObjectParent(ts,cvode->ydot);CHKERRQ(ierr); 3187cb94ee6SHong Zhang 3197cb94ee6SHong Zhang /* 3207cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 3217cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 3227cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 3237cb94ee6SHong Zhang */ 324778a2246SBarry Smith ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,1,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr); 325778a2246SBarry Smith ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,1,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr); 3267cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr); 3277cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr); 32816016685SHong Zhang 32916016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 33016016685SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 331e32f2f54SBarry Smith if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails"); 33216016685SHong Zhang cvode->mem = mem; 33316016685SHong Zhang 33416016685SHong Zhang /* Set the pointer to user-defined data */ 33516016685SHong Zhang flag = CVodeSetUserData(mem, ts); 336e32f2f54SBarry Smith if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails"); 33716016685SHong Zhang 338fc6b9e64SJed Brown /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 339cdbf8f93SLisandro Dalcin flag = CVodeSetInitStep(mem,(realtype)ts->time_step); 340fc6b9e64SJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed"); 341f1cd61daSJed Brown if (cvode->mindt > 0) { 342f1cd61daSJed Brown flag = CVodeSetMinStep(mem,(realtype)cvode->mindt); 3439f94935aSHong Zhang if (flag){ 3449f94935aSHong Zhang if (flag == CV_MEM_NULL){ 3459f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 3469f94935aSHong Zhang } else if (flag == CV_ILL_INPUT){ 3479f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size"); 3489f94935aSHong Zhang } else { 3499f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed"); 3509f94935aSHong Zhang } 3519f94935aSHong Zhang } 352f1cd61daSJed Brown } 353f1cd61daSJed Brown if (cvode->maxdt > 0) { 354f1cd61daSJed Brown flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt); 355f1cd61daSJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed"); 356f1cd61daSJed Brown } 357f1cd61daSJed Brown 35816016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 35916016685SHong Zhang * user's right hand side function in u'=f(t,u), the inital time T0, and 36016016685SHong Zhang * the initial dependent variable vector cvode->y */ 36116016685SHong Zhang flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 36216016685SHong Zhang if (flag){ 363e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag); 36416016685SHong Zhang } 36516016685SHong Zhang 3669f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 36716016685SHong Zhang flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 36816016685SHong Zhang if (flag){ 369e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag); 37016016685SHong Zhang } 37116016685SHong Zhang 3729f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 3739f94935aSHong Zhang flag = CVodeSetMaxNumSteps(mem,ts->max_steps); 37416016685SHong Zhang 37516016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 37616016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 377dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr); 378dcbc6d53SSean Farley ierr = PCGetType(pc,&pctype);CHKERRQ(ierr); 379251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr); 38016016685SHong Zhang if (pcnone){ 38116016685SHong Zhang flag = CVSpgmr(mem,PREC_NONE,0); 382e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 38316016685SHong Zhang } else { 384f61b2b6aSHong Zhang flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl); 385e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 38616016685SHong Zhang 38716016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 38816016685SHong Zhang and the pointer to the user-defined block data */ 38916016685SHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 390e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 39116016685SHong Zhang } 39216016685SHong Zhang 39316016685SHong Zhang flag = CVSpilsSetGSType(mem, MODIFIED_GS); 394e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag); 3957cb94ee6SHong Zhang PetscFunctionReturn(0); 3967cb94ee6SHong Zhang } 3977cb94ee6SHong Zhang 3986fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 3996a6fc655SJed Brown const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",0}; 4006fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 4016a6fc655SJed Brown const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",0}; 402a04cf4d8SBarry Smith 4037cb94ee6SHong Zhang #undef __FUNCT__ 404214bc6a2SJed Brown #define __FUNCT__ "TSSetFromOptions_Sundials" 405214bc6a2SJed Brown PetscErrorCode TSSetFromOptions_Sundials(TS ts) 4067cb94ee6SHong Zhang { 4077cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4087cb94ee6SHong Zhang PetscErrorCode ierr; 4097cb94ee6SHong Zhang int indx; 410ace3abfcSBarry Smith PetscBool flag; 4117bda3a07SJed Brown PC pc; 4127cb94ee6SHong Zhang 4137cb94ee6SHong Zhang PetscFunctionBegin; 4147cb94ee6SHong Zhang ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr); 4156fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr); 4167cb94ee6SHong Zhang if (flag) { 4176fadb2cdSHong Zhang ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr); 4187cb94ee6SHong Zhang } 4196fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr); 4207cb94ee6SHong Zhang if (flag) { 4217cb94ee6SHong Zhang ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 4227cb94ee6SHong Zhang } 4237cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr); 4247cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr); 425f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr); 426f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr); 4277cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr); 428f61b2b6aSHong Zhang ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,&flag);CHKERRQ(ierr); 429acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr); 4307cb94ee6SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4317bda3a07SJed Brown ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 4327bda3a07SJed Brown ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 4337cb94ee6SHong Zhang PetscFunctionReturn(0); 4347cb94ee6SHong Zhang } 4357cb94ee6SHong Zhang 4367cb94ee6SHong Zhang #undef __FUNCT__ 4377cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials" 4387cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 4397cb94ee6SHong Zhang { 4407cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4417cb94ee6SHong Zhang PetscErrorCode ierr; 4427cb94ee6SHong Zhang char *type; 4437cb94ee6SHong Zhang char atype[] = "Adams"; 4447cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 445ace3abfcSBarry Smith PetscBool iascii,isstring; 4462c823083SHong Zhang long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 4472c823083SHong Zhang PetscInt qlast,qcur; 4485d47aee6SHong Zhang PetscReal hinused,hlast,hcur,tcur,tolsfac; 44942528757SHong Zhang PC pc; 4507cb94ee6SHong Zhang 4517cb94ee6SHong Zhang PetscFunctionBegin; 4527cb94ee6SHong Zhang if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;} 4537cb94ee6SHong Zhang else {type = btype;} 4547cb94ee6SHong Zhang 455251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 456251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 4577cb94ee6SHong Zhang if (iascii) { 4587cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr); 4597cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr); 4607cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr); 4617cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr); 462f61b2b6aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr); 4637cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 4647cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4657cb94ee6SHong Zhang } else { 4667cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4677cb94ee6SHong Zhang } 468f1cd61daSJed Brown if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);} 469f1cd61daSJed Brown if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);} 4702c823083SHong Zhang 4715d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4725d47aee6SHong Zhang ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr); 4735d47aee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr); 4742c823083SHong Zhang ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals, 4752c823083SHong Zhang &nlinsetups,&nfails,&qlast,&qcur, 4762c823083SHong Zhang &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr); 4772c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr); 4782c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr); 4792c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr); 4802c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr); 4812c823083SHong Zhang 4822c823083SHong Zhang ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr); 4832c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr); 4842c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr); 4852c823083SHong Zhang 4862c823083SHong Zhang ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */ 4872c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr); 4882c823083SHong Zhang ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr); 4892c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr); 49042528757SHong Zhang 49142528757SHong Zhang ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr); 49242528757SHong Zhang ierr = PCView(pc,viewer);CHKERRQ(ierr); 4932c823083SHong Zhang ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4942c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr); 4952c823083SHong Zhang ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr); 4962c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr); 49742528757SHong Zhang 4982c823083SHong Zhang ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4992c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr); 5002c823083SHong Zhang ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr); 5012c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr); 5027cb94ee6SHong Zhang } else if (isstring) { 5037cb94ee6SHong Zhang ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 5047cb94ee6SHong Zhang } 5057cb94ee6SHong Zhang PetscFunctionReturn(0); 5067cb94ee6SHong Zhang } 5077cb94ee6SHong Zhang 5087cb94ee6SHong Zhang 5097cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 5107cb94ee6SHong Zhang EXTERN_C_BEGIN 5117cb94ee6SHong Zhang #undef __FUNCT__ 5127cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials" 5137087cfbeSBarry Smith PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 5147cb94ee6SHong Zhang { 5157cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5167cb94ee6SHong Zhang 5177cb94ee6SHong Zhang PetscFunctionBegin; 5187cb94ee6SHong Zhang cvode->cvode_type = type; 5197cb94ee6SHong Zhang PetscFunctionReturn(0); 5207cb94ee6SHong Zhang } 5217cb94ee6SHong Zhang EXTERN_C_END 5227cb94ee6SHong Zhang 5237cb94ee6SHong Zhang EXTERN_C_BEGIN 5247cb94ee6SHong Zhang #undef __FUNCT__ 525f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl_Sundials" 526f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl) 5277cb94ee6SHong Zhang { 5287cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5297cb94ee6SHong Zhang 5307cb94ee6SHong Zhang PetscFunctionBegin; 531f61b2b6aSHong Zhang cvode->maxl = maxl; 5327cb94ee6SHong Zhang PetscFunctionReturn(0); 5337cb94ee6SHong Zhang } 5347cb94ee6SHong Zhang EXTERN_C_END 5357cb94ee6SHong Zhang 5367cb94ee6SHong Zhang EXTERN_C_BEGIN 5377cb94ee6SHong Zhang #undef __FUNCT__ 5387cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials" 5397087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,double tol) 5407cb94ee6SHong Zhang { 5417cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5427cb94ee6SHong Zhang 5437cb94ee6SHong Zhang PetscFunctionBegin; 5447cb94ee6SHong Zhang cvode->linear_tol = tol; 5457cb94ee6SHong Zhang PetscFunctionReturn(0); 5467cb94ee6SHong Zhang } 5477cb94ee6SHong Zhang EXTERN_C_END 5487cb94ee6SHong Zhang 5497cb94ee6SHong Zhang EXTERN_C_BEGIN 5507cb94ee6SHong Zhang #undef __FUNCT__ 5517cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials" 5527087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 5537cb94ee6SHong Zhang { 5547cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5557cb94ee6SHong Zhang 5567cb94ee6SHong Zhang PetscFunctionBegin; 5577cb94ee6SHong Zhang cvode->gtype = type; 5587cb94ee6SHong Zhang PetscFunctionReturn(0); 5597cb94ee6SHong Zhang } 5607cb94ee6SHong Zhang EXTERN_C_END 5617cb94ee6SHong Zhang 5627cb94ee6SHong Zhang EXTERN_C_BEGIN 5637cb94ee6SHong Zhang #undef __FUNCT__ 5647cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials" 5657087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel) 5667cb94ee6SHong Zhang { 5677cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5687cb94ee6SHong Zhang 5697cb94ee6SHong Zhang PetscFunctionBegin; 5707cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 5717cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 5727cb94ee6SHong Zhang PetscFunctionReturn(0); 5737cb94ee6SHong Zhang } 5747cb94ee6SHong Zhang EXTERN_C_END 5757cb94ee6SHong Zhang 5767cb94ee6SHong Zhang EXTERN_C_BEGIN 5777cb94ee6SHong Zhang #undef __FUNCT__ 578f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials" 5797087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 580f1cd61daSJed Brown { 581f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 582f1cd61daSJed Brown 583f1cd61daSJed Brown PetscFunctionBegin; 584f1cd61daSJed Brown cvode->mindt = mindt; 585f1cd61daSJed Brown PetscFunctionReturn(0); 586f1cd61daSJed Brown } 587f1cd61daSJed Brown EXTERN_C_END 588f1cd61daSJed Brown 589f1cd61daSJed Brown EXTERN_C_BEGIN 590f1cd61daSJed Brown #undef __FUNCT__ 591f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials" 5927087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 593f1cd61daSJed Brown { 594f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 595f1cd61daSJed Brown 596f1cd61daSJed Brown PetscFunctionBegin; 597f1cd61daSJed Brown cvode->maxdt = maxdt; 598f1cd61daSJed Brown PetscFunctionReturn(0); 599f1cd61daSJed Brown } 600f1cd61daSJed Brown EXTERN_C_END 601f1cd61daSJed Brown 602f1cd61daSJed Brown EXTERN_C_BEGIN 603f1cd61daSJed Brown #undef __FUNCT__ 6047cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials" 6057087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc) 6067cb94ee6SHong Zhang { 607dcbc6d53SSean Farley SNES snes; 608dcbc6d53SSean Farley KSP ksp; 609dcbc6d53SSean Farley PetscErrorCode ierr; 6107cb94ee6SHong Zhang 6117cb94ee6SHong Zhang PetscFunctionBegin; 612dcbc6d53SSean Farley ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 613dcbc6d53SSean Farley ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 614dcbc6d53SSean Farley ierr = KSPGetPC(ksp,pc); CHKERRQ(ierr); 6157cb94ee6SHong Zhang PetscFunctionReturn(0); 6167cb94ee6SHong Zhang } 6177cb94ee6SHong Zhang EXTERN_C_END 6187cb94ee6SHong Zhang 6197cb94ee6SHong Zhang EXTERN_C_BEGIN 6207cb94ee6SHong Zhang #undef __FUNCT__ 6217cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials" 6227087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 6237cb94ee6SHong Zhang { 6247cb94ee6SHong Zhang PetscFunctionBegin; 6255ef26d82SJed Brown if (nonlin) *nonlin = ts->snes_its; 6265ef26d82SJed Brown if (lin) *lin = ts->ksp_its; 6277cb94ee6SHong Zhang PetscFunctionReturn(0); 6287cb94ee6SHong Zhang } 6297cb94ee6SHong Zhang EXTERN_C_END 6307cb94ee6SHong Zhang 6317cb94ee6SHong Zhang EXTERN_C_BEGIN 6327cb94ee6SHong Zhang #undef __FUNCT__ 6332bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials" 6347087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s) 6352bfc04deSHong Zhang { 6362bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 6372bfc04deSHong Zhang 6382bfc04deSHong Zhang PetscFunctionBegin; 6392bfc04deSHong Zhang cvode->monitorstep = s; 6402bfc04deSHong Zhang PetscFunctionReturn(0); 6412bfc04deSHong Zhang } 6422bfc04deSHong Zhang EXTERN_C_END 6437cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 6447cb94ee6SHong Zhang 6457cb94ee6SHong Zhang #undef __FUNCT__ 6467cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations" 6477cb94ee6SHong Zhang /*@C 6487cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 6497cb94ee6SHong Zhang 6507cb94ee6SHong Zhang Not Collective 6517cb94ee6SHong Zhang 6527cb94ee6SHong Zhang Input parameters: 6537cb94ee6SHong Zhang . ts - the time-step context 6547cb94ee6SHong Zhang 6557cb94ee6SHong Zhang Output Parameters: 6567cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 6577cb94ee6SHong Zhang - lin - number of linear iterations 6587cb94ee6SHong Zhang 6597cb94ee6SHong Zhang Level: advanced 6607cb94ee6SHong Zhang 6617cb94ee6SHong Zhang Notes: 6627cb94ee6SHong Zhang These return the number since the creation of the TS object 6637cb94ee6SHong Zhang 6647cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations 6657cb94ee6SHong Zhang 666f61b2b6aSHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetMaxl(), 6677cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 668f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 669a43b19c4SJed Brown TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime() 6707cb94ee6SHong Zhang 6717cb94ee6SHong Zhang @*/ 6727087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 6737cb94ee6SHong Zhang { 6744ac538c5SBarry Smith PetscErrorCode ierr; 6757cb94ee6SHong Zhang 6767cb94ee6SHong Zhang PetscFunctionBegin; 6774ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr); 6787cb94ee6SHong Zhang PetscFunctionReturn(0); 6797cb94ee6SHong Zhang } 6807cb94ee6SHong Zhang 6817cb94ee6SHong Zhang #undef __FUNCT__ 6827cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType" 6837cb94ee6SHong Zhang /*@ 6847cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 6857cb94ee6SHong Zhang 6863f9fe445SBarry Smith Logically Collective on TS 6877cb94ee6SHong Zhang 6887cb94ee6SHong Zhang Input parameters: 6897cb94ee6SHong Zhang + ts - the time-step context 6907cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 6917cb94ee6SHong Zhang 6927cb94ee6SHong Zhang Level: intermediate 6937cb94ee6SHong Zhang 6947cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula 6957cb94ee6SHong Zhang 696f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetMaxl(), 6977cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 698f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 6997cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 700a43b19c4SJed Brown TSSetExactFinalTime() 7017cb94ee6SHong Zhang @*/ 7027087cfbeSBarry Smith PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type) 7037cb94ee6SHong Zhang { 7044ac538c5SBarry Smith PetscErrorCode ierr; 7057cb94ee6SHong Zhang 7067cb94ee6SHong Zhang PetscFunctionBegin; 7074ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr); 7087cb94ee6SHong Zhang PetscFunctionReturn(0); 7097cb94ee6SHong Zhang } 7107cb94ee6SHong Zhang 7117cb94ee6SHong Zhang #undef __FUNCT__ 712f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl" 7137cb94ee6SHong Zhang /*@ 714f61b2b6aSHong Zhang TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 7157cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 716f61b2b6aSHong Zhang this is the maximum number of GMRES steps that will be used. 7177cb94ee6SHong Zhang 7183f9fe445SBarry Smith Logically Collective on TS 7197cb94ee6SHong Zhang 7207cb94ee6SHong Zhang Input parameters: 7217cb94ee6SHong Zhang + ts - the time-step context 722f61b2b6aSHong Zhang - maxl - number of direction vectors (the dimension of Krylov subspace). 7237cb94ee6SHong Zhang 7247cb94ee6SHong Zhang Level: advanced 7257cb94ee6SHong Zhang 726f61b2b6aSHong Zhang .keywords: GMRES 7277cb94ee6SHong Zhang 7287cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 7297cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 730f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7317cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 732a43b19c4SJed Brown TSSetExactFinalTime() 7337cb94ee6SHong Zhang 7347cb94ee6SHong Zhang @*/ 735f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl(TS ts,PetscInt maxl) 7367cb94ee6SHong Zhang { 7374ac538c5SBarry Smith PetscErrorCode ierr; 7387cb94ee6SHong Zhang 7397cb94ee6SHong Zhang PetscFunctionBegin; 740f61b2b6aSHong Zhang PetscValidLogicalCollectiveInt(ts,maxl,2); 741f61b2b6aSHong Zhang ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));CHKERRQ(ierr); 7427cb94ee6SHong Zhang PetscFunctionReturn(0); 7437cb94ee6SHong Zhang } 7447cb94ee6SHong Zhang 7457cb94ee6SHong Zhang #undef __FUNCT__ 7467cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance" 7477cb94ee6SHong Zhang /*@ 7487cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 7497cb94ee6SHong Zhang system by SUNDIALS. 7507cb94ee6SHong Zhang 7513f9fe445SBarry Smith Logically Collective on TS 7527cb94ee6SHong Zhang 7537cb94ee6SHong Zhang Input parameters: 7547cb94ee6SHong Zhang + ts - the time-step context 7557cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 7567cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 7577cb94ee6SHong Zhang 7587cb94ee6SHong Zhang Level: advanced 7597cb94ee6SHong Zhang 7607cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS 7617cb94ee6SHong Zhang 762f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7637cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 764f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7657cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 766a43b19c4SJed Brown TSSetExactFinalTime() 7677cb94ee6SHong Zhang 7687cb94ee6SHong Zhang @*/ 7697087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance(TS ts,double tol) 7707cb94ee6SHong Zhang { 7714ac538c5SBarry Smith PetscErrorCode ierr; 7727cb94ee6SHong Zhang 7737cb94ee6SHong Zhang PetscFunctionBegin; 774c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,tol,2); 7754ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr); 7767cb94ee6SHong Zhang PetscFunctionReturn(0); 7777cb94ee6SHong Zhang } 7787cb94ee6SHong Zhang 7797cb94ee6SHong Zhang #undef __FUNCT__ 7807cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType" 7817cb94ee6SHong Zhang /*@ 7827cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 7837cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 7847cb94ee6SHong Zhang 7853f9fe445SBarry Smith Logically Collective on TS 7867cb94ee6SHong Zhang 7877cb94ee6SHong Zhang Input parameters: 7887cb94ee6SHong Zhang + ts - the time-step context 7897cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 7907cb94ee6SHong Zhang 7917cb94ee6SHong Zhang Level: advanced 7927cb94ee6SHong Zhang 7937cb94ee6SHong Zhang .keywords: Sundials, orthogonalization 7947cb94ee6SHong Zhang 795f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7967cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 797f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7987cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 799a43b19c4SJed Brown TSSetExactFinalTime() 8007cb94ee6SHong Zhang 8017cb94ee6SHong Zhang @*/ 8027087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 8037cb94ee6SHong Zhang { 8044ac538c5SBarry Smith PetscErrorCode ierr; 8057cb94ee6SHong Zhang 8067cb94ee6SHong Zhang PetscFunctionBegin; 8074ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr); 8087cb94ee6SHong Zhang PetscFunctionReturn(0); 8097cb94ee6SHong Zhang } 8107cb94ee6SHong Zhang 8117cb94ee6SHong Zhang #undef __FUNCT__ 8127cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance" 8137cb94ee6SHong Zhang /*@ 8147cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 8157cb94ee6SHong Zhang Sundials for error control. 8167cb94ee6SHong Zhang 8173f9fe445SBarry Smith Logically Collective on TS 8187cb94ee6SHong Zhang 8197cb94ee6SHong Zhang Input parameters: 8207cb94ee6SHong Zhang + ts - the time-step context 8217cb94ee6SHong Zhang . aabs - the absolute tolerance 8227cb94ee6SHong Zhang - rel - the relative tolerance 8237cb94ee6SHong Zhang 8247cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 8257cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 8267cb94ee6SHong Zhang 8277cb94ee6SHong Zhang Level: intermediate 8287cb94ee6SHong Zhang 8297cb94ee6SHong Zhang .keywords: Sundials, tolerance 8307cb94ee6SHong Zhang 831f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(), 8327cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 833f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8347cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 835a43b19c4SJed Brown TSSetExactFinalTime() 8367cb94ee6SHong Zhang 8377cb94ee6SHong Zhang @*/ 8387087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance(TS ts,double aabs,double rel) 8397cb94ee6SHong Zhang { 8404ac538c5SBarry Smith PetscErrorCode ierr; 8417cb94ee6SHong Zhang 8427cb94ee6SHong Zhang PetscFunctionBegin; 8434ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr); 8447cb94ee6SHong Zhang PetscFunctionReturn(0); 8457cb94ee6SHong Zhang } 8467cb94ee6SHong Zhang 8477cb94ee6SHong Zhang #undef __FUNCT__ 8487cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC" 8497cb94ee6SHong Zhang /*@ 8507cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 8517cb94ee6SHong Zhang 8527cb94ee6SHong Zhang Input Parameter: 8537cb94ee6SHong Zhang . ts - the time-step context 8547cb94ee6SHong Zhang 8557cb94ee6SHong Zhang Output Parameter: 8567cb94ee6SHong Zhang . pc - the preconditioner context 8577cb94ee6SHong Zhang 8587cb94ee6SHong Zhang Level: advanced 8597cb94ee6SHong Zhang 860f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 8617cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 862f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8637cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 8647cb94ee6SHong Zhang @*/ 8657087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC(TS ts,PC *pc) 8667cb94ee6SHong Zhang { 8674ac538c5SBarry Smith PetscErrorCode ierr; 8687cb94ee6SHong Zhang 8697cb94ee6SHong Zhang PetscFunctionBegin; 8704ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr); 8717cb94ee6SHong Zhang PetscFunctionReturn(0); 8727cb94ee6SHong Zhang } 8737cb94ee6SHong Zhang 8747cb94ee6SHong Zhang #undef __FUNCT__ 875f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep" 876f1cd61daSJed Brown /*@ 877f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 878f1cd61daSJed Brown 879f1cd61daSJed Brown Input Parameter: 880f1cd61daSJed Brown + ts - the time-step context 881f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 882f1cd61daSJed Brown 883fc6b9e64SJed Brown Note: 884fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 885fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 886fc6b9e64SJed Brown 887f1cd61daSJed Brown Level: beginner 888f1cd61daSJed Brown 889f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 890f1cd61daSJed Brown @*/ 8917087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 892f1cd61daSJed Brown { 8934ac538c5SBarry Smith PetscErrorCode ierr; 894f1cd61daSJed Brown 895f1cd61daSJed Brown PetscFunctionBegin; 8964ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr); 897f1cd61daSJed Brown PetscFunctionReturn(0); 898f1cd61daSJed Brown } 899f1cd61daSJed Brown 900f1cd61daSJed Brown #undef __FUNCT__ 901f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep" 902f1cd61daSJed Brown /*@ 903f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 904f1cd61daSJed Brown 905f1cd61daSJed Brown Input Parameter: 906f1cd61daSJed Brown + ts - the time-step context 907f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 908f1cd61daSJed Brown 909f1cd61daSJed Brown Level: beginner 910f1cd61daSJed Brown 911f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 912f1cd61daSJed Brown @*/ 9137087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 914f1cd61daSJed Brown { 9154ac538c5SBarry Smith PetscErrorCode ierr; 916f1cd61daSJed Brown 917f1cd61daSJed Brown PetscFunctionBegin; 9184ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 919f1cd61daSJed Brown PetscFunctionReturn(0); 920f1cd61daSJed Brown } 921f1cd61daSJed Brown 922f1cd61daSJed Brown #undef __FUNCT__ 9232bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps" 9242bfc04deSHong Zhang /*@ 9252bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 9262bfc04deSHong Zhang 9272bfc04deSHong Zhang Input Parameter: 9282bfc04deSHong Zhang + ts - the time-step context 9292bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 9302bfc04deSHong Zhang 9312bfc04deSHong Zhang Level: beginner 9322bfc04deSHong Zhang 933f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 9342bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 935f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 9362bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 9372bfc04deSHong Zhang @*/ 9387087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft) 9392bfc04deSHong Zhang { 9404ac538c5SBarry Smith PetscErrorCode ierr; 9412bfc04deSHong Zhang 9422bfc04deSHong Zhang PetscFunctionBegin; 9434ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 9442bfc04deSHong Zhang PetscFunctionReturn(0); 9452bfc04deSHong Zhang } 9467cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 9477cb94ee6SHong Zhang /*MC 94896f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 9497cb94ee6SHong Zhang 9507cb94ee6SHong Zhang Options Database: 9517cb94ee6SHong Zhang + -ts_sundials_type <bdf,adams> 9527cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 9537cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 9547cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 9557cb94ee6SHong Zhang . -ts_sundials_linear_tolerance <tol> 956f61b2b6aSHong Zhang . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 95716016685SHong Zhang - -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps 95816016685SHong Zhang 9597cb94ee6SHong Zhang 9607cb94ee6SHong Zhang Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply 9617cb94ee6SHong Zhang only PETSc PC options 9627cb94ee6SHong Zhang 9637cb94ee6SHong Zhang Level: beginner 9647cb94ee6SHong Zhang 965f61b2b6aSHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(), 966a43b19c4SJed Brown TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime() 9677cb94ee6SHong Zhang 9687cb94ee6SHong Zhang M*/ 9697cb94ee6SHong Zhang EXTERN_C_BEGIN 9707cb94ee6SHong Zhang #undef __FUNCT__ 9717cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials" 9727087cfbeSBarry Smith PetscErrorCode TSCreate_Sundials(TS ts) 9737cb94ee6SHong Zhang { 9747cb94ee6SHong Zhang TS_Sundials *cvode; 9757cb94ee6SHong Zhang PetscErrorCode ierr; 97642528757SHong Zhang PC pc; 9777cb94ee6SHong Zhang 9787cb94ee6SHong Zhang PetscFunctionBegin; 979277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Sundials; 98028adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 98128adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 982214bc6a2SJed Brown ts->ops->setup = TSSetUp_Sundials; 983b4eba00bSSean Farley ts->ops->step = TSStep_Sundials; 984b4eba00bSSean Farley ts->ops->interpolate = TSInterpolate_Sundials; 985214bc6a2SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Sundials; 9867cb94ee6SHong Zhang 98738f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr); 9887cb94ee6SHong Zhang ts->data = (void*)cvode; 9896fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 9906fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 991f61b2b6aSHong Zhang cvode->maxl = 5; 9927cb94ee6SHong Zhang cvode->linear_tol = .05; 9937cb94ee6SHong Zhang 994b4eba00bSSean Farley cvode->monitorstep = PETSC_TRUE; 9957cb94ee6SHong Zhang 996a07356b8SSatish Balay ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr); 997f1cd61daSJed Brown 998f1cd61daSJed Brown cvode->mindt = -1.; 999f1cd61daSJed Brown cvode->maxdt = -1.; 1000f1cd61daSJed Brown 10017cb94ee6SHong Zhang /* set tolerance for Sundials */ 10027cb94ee6SHong Zhang cvode->reltol = 1e-6; 10032c823083SHong Zhang cvode->abstol = 1e-6; 10047cb94ee6SHong Zhang 100542528757SHong Zhang /* set PCNONE as default pctype */ 100642528757SHong Zhang ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr); 100742528757SHong Zhang ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 100842528757SHong Zhang 1009b4eba00bSSean Farley if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE; 1010a43b19c4SJed Brown 10117cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials", 10127cb94ee6SHong Zhang TSSundialsSetType_Sundials);CHKERRQ(ierr); 1013f61b2b6aSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C", 1014f61b2b6aSHong Zhang "TSSundialsSetMaxl_Sundials", 1015f61b2b6aSHong Zhang TSSundialsSetMaxl_Sundials);CHKERRQ(ierr); 10167cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C", 10177cb94ee6SHong Zhang "TSSundialsSetLinearTolerance_Sundials", 10187cb94ee6SHong Zhang TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 10197cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C", 10207cb94ee6SHong Zhang "TSSundialsSetGramSchmidtType_Sundials", 10217cb94ee6SHong Zhang TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 10227cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C", 10237cb94ee6SHong Zhang "TSSundialsSetTolerance_Sundials", 10247cb94ee6SHong Zhang TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 1025f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C", 1026f1cd61daSJed Brown "TSSundialsSetMinTimeStep_Sundials", 1027f1cd61daSJed Brown TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr); 1028f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C", 1029f1cd61daSJed Brown "TSSundialsSetMaxTimeStep_Sundials", 1030f1cd61daSJed Brown TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr); 10317cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C", 10327cb94ee6SHong Zhang "TSSundialsGetPC_Sundials", 10337cb94ee6SHong Zhang TSSundialsGetPC_Sundials);CHKERRQ(ierr); 10347cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C", 10357cb94ee6SHong Zhang "TSSundialsGetIterations_Sundials", 10367cb94ee6SHong Zhang TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 10372bfc04deSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C", 10382bfc04deSHong Zhang "TSSundialsMonitorInternalSteps_Sundials", 10392bfc04deSHong Zhang TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr); 10402bfc04deSHong Zhang 10417cb94ee6SHong Zhang PetscFunctionReturn(0); 10427cb94ee6SHong Zhang } 10437cb94ee6SHong Zhang EXTERN_C_END 1044