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; 847adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)ts)->comm; 857cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 860679a0aeSJed 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); 964ac7836bSHong Zhang 977cb94ee6SHong Zhang /* now compute the right hand side function */ 98*bc0cc02bSJed Brown if (!ts->userops->ifunction) { 99*bc0cc02bSJed Brown ierr = TSComputeRHSFunction(ts,t,yy,yyd);CHKERRQ(ierr); 100*bc0cc02bSJed Brown } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */ 101*bc0cc02bSJed Brown ierr = VecZeroEntries(yydot);CHKERRQ(ierr); 1020679a0aeSJed Brown ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE); CHKERRABORT(comm,ierr); 1030679a0aeSJed Brown ierr = VecScale(yyd,-1.);CHKERRQ(ierr); 104*bc0cc02bSJed Brown } 1057cb94ee6SHong Zhang ierr = VecResetArray(yy); CHKERRABORT(comm,ierr); 1067cb94ee6SHong Zhang ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr); 1074ac7836bSHong Zhang PetscFunctionReturn(0); 1087cb94ee6SHong Zhang } 1097cb94ee6SHong Zhang 1107cb94ee6SHong Zhang /* 111b4eba00bSSean Farley TSStep_Sundials - Calls Sundials to integrate the ODE. 1127cb94ee6SHong Zhang */ 1137cb94ee6SHong Zhang #undef __FUNCT__ 114b4eba00bSSean Farley #define __FUNCT__ "TSStep_Sundials" 115b4eba00bSSean Farley PetscErrorCode TSStep_Sundials(TS ts) 1167cb94ee6SHong Zhang { 1177cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1187cb94ee6SHong Zhang PetscErrorCode ierr; 119b4eba00bSSean Farley PetscInt flag; 1209f94935aSHong Zhang long int its,nsteps; 1217cb94ee6SHong Zhang realtype t,tout; 1227cb94ee6SHong Zhang PetscScalar *y_data; 1237cb94ee6SHong Zhang void *mem; 1247cb94ee6SHong Zhang 1257cb94ee6SHong Zhang PetscFunctionBegin; 12616016685SHong Zhang mem = cvode->mem; 1277cb94ee6SHong Zhang tout = ts->max_time; 1287cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr); 1297cb94ee6SHong Zhang N_VSetArrayPointer((realtype *)y_data,cvode->y); 1307cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 131186e87acSLisandro Dalcin 1323f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 133186e87acSLisandro Dalcin 1342bfc04deSHong Zhang if (cvode->monitorstep) { 13528adb3f7SHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP); 1362bfc04deSHong Zhang } else { 1372bfc04deSHong Zhang flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL); 1382bfc04deSHong Zhang } 1399f94935aSHong Zhang 1409f94935aSHong Zhang if (flag){ /* display error message */ 1419f94935aSHong Zhang switch (flag){ 1429f94935aSHong Zhang case CV_ILL_INPUT: 1439f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT"); 1449f94935aSHong Zhang break; 1459f94935aSHong Zhang case CV_TOO_CLOSE: 1469f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE"); 1479f94935aSHong Zhang break; 1483c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1499f94935aSHong Zhang PetscReal tcur; 1509f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 1519f94935aSHong Zhang ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr); 1529f94935aSHong 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); 1533c7fefeeSJed Brown } break; 1549f94935aSHong Zhang case CV_TOO_MUCH_ACC: 1559f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC"); 1569f94935aSHong Zhang break; 1579f94935aSHong Zhang case CV_ERR_FAILURE: 1589f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE"); 1599f94935aSHong Zhang break; 1609f94935aSHong Zhang case CV_CONV_FAILURE: 1619f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE"); 1629f94935aSHong Zhang break; 1639f94935aSHong Zhang case CV_LINIT_FAIL: 1649f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL"); 1659f94935aSHong Zhang break; 1669f94935aSHong Zhang case CV_LSETUP_FAIL: 1679f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL"); 1689f94935aSHong Zhang break; 1699f94935aSHong Zhang case CV_LSOLVE_FAIL: 1709f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL"); 1719f94935aSHong Zhang break; 1729f94935aSHong Zhang case CV_RHSFUNC_FAIL: 1739f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL"); 1749f94935aSHong Zhang break; 1759f94935aSHong Zhang case CV_FIRST_RHSFUNC_ERR: 1769f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR"); 1779f94935aSHong Zhang break; 1789f94935aSHong Zhang case CV_REPTD_RHSFUNC_ERR: 1799f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR"); 1809f94935aSHong Zhang break; 1819f94935aSHong Zhang case CV_UNREC_RHSFUNC_ERR: 1829f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR"); 1839f94935aSHong Zhang break; 1849f94935aSHong Zhang case CV_RTFUNC_FAIL: 1859f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL"); 1869f94935aSHong Zhang break; 1879f94935aSHong Zhang default: 1889f94935aSHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag); 1899f94935aSHong Zhang } 1909f94935aSHong Zhang } 1919f94935aSHong Zhang 1927cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 1937cb94ee6SHong Zhang ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr); 1947cb94ee6SHong Zhang ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr); 1957cb94ee6SHong Zhang ierr = VecResetArray(cvode->w1); CHKERRQ(ierr); 196*bc0cc02bSJed Brown ierr = VecCopy(cvode->update,ts->vec_sol);CHKERRQ(ierr); 1977cb94ee6SHong Zhang ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr); 1984ac7836bSHong Zhang ierr = CVSpilsGetNumLinIters(mem, &its); 199186e87acSLisandro Dalcin ts->nonlinear_its = its; ts->linear_its = its; 200186e87acSLisandro Dalcin 201186e87acSLisandro Dalcin ts->time_step = t - ts->ptime; 202186e87acSLisandro Dalcin ts->ptime = t; 2037cb94ee6SHong Zhang ts->steps++; 204186e87acSLisandro Dalcin 2059f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 206b4eba00bSSean Farley if (!cvode->monitorstep) ts->steps = nsteps; 207b4eba00bSSean Farley PetscFunctionReturn(0); 208b4eba00bSSean Farley } 209b4eba00bSSean Farley 210b4eba00bSSean Farley #undef __FUNCT__ 211b4eba00bSSean Farley #define __FUNCT__ "TSInterpolate_Sundials" 212b4eba00bSSean Farley static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X) 213b4eba00bSSean Farley { 214b4eba00bSSean Farley TS_Sundials *cvode = (TS_Sundials*)ts->data; 215b4eba00bSSean Farley N_Vector y; 216b4eba00bSSean Farley PetscErrorCode ierr; 217b4eba00bSSean Farley PetscScalar *x_data; 218b4eba00bSSean Farley PetscInt glosize,locsize; 219b4eba00bSSean Farley 220b4eba00bSSean Farley PetscFunctionBegin; 221b4eba00bSSean Farley 222b4eba00bSSean Farley /* get the vector size */ 223b4eba00bSSean Farley ierr = VecGetSize(X,&glosize);CHKERRQ(ierr); 224b4eba00bSSean Farley ierr = VecGetLocalSize(X,&locsize);CHKERRQ(ierr); 225b4eba00bSSean Farley 226b4eba00bSSean Farley /* allocate the memory for N_Vec y */ 227b4eba00bSSean Farley y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 228b4eba00bSSean Farley if (!y) SETERRQ(PETSC_COMM_SELF,1,"Interpolated y is not allocated"); 229b4eba00bSSean Farley 230b4eba00bSSean Farley ierr = VecGetArray(X,&x_data);CHKERRQ(ierr); 231b4eba00bSSean Farley N_VSetArrayPointer((realtype *)x_data,y); 232b4eba00bSSean Farley ierr = CVodeGetDky(cvode->mem,t,0,y);CHKERRQ(ierr); 233b4eba00bSSean Farley ierr = VecRestoreArray(X,&x_data);CHKERRQ(ierr); 234b4eba00bSSean Farley 2357cb94ee6SHong Zhang PetscFunctionReturn(0); 2367cb94ee6SHong Zhang } 2377cb94ee6SHong Zhang 2387cb94ee6SHong Zhang #undef __FUNCT__ 239277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Sundials" 240277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts) 241277b19d0SLisandro Dalcin { 242277b19d0SLisandro Dalcin TS_Sundials *cvode = (TS_Sundials*)ts->data; 243277b19d0SLisandro Dalcin PetscErrorCode ierr; 244277b19d0SLisandro Dalcin 245277b19d0SLisandro Dalcin PetscFunctionBegin; 2465c6b4a3dSJed Brown ierr = VecDestroy(&cvode->update);CHKERRQ(ierr); 2470679a0aeSJed Brown ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr); 2485c6b4a3dSJed Brown ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr); 2495c6b4a3dSJed Brown ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr); 250277b19d0SLisandro Dalcin if (cvode->mem) {CVodeFree(&cvode->mem);} 251277b19d0SLisandro Dalcin PetscFunctionReturn(0); 252277b19d0SLisandro Dalcin } 253277b19d0SLisandro Dalcin 254277b19d0SLisandro Dalcin #undef __FUNCT__ 2557cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials" 2567cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts) 2577cb94ee6SHong Zhang { 2587cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2597cb94ee6SHong Zhang PetscErrorCode ierr; 2607cb94ee6SHong Zhang 2617cb94ee6SHong Zhang PetscFunctionBegin; 262277b19d0SLisandro Dalcin ierr = TSReset_Sundials(ts);CHKERRQ(ierr); 263a07356b8SSatish Balay ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr); 264277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 265335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","",PETSC_NULL);CHKERRQ(ierr); 266f61b2b6aSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C","",PETSC_NULL);CHKERRQ(ierr); 267335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",PETSC_NULL);CHKERRQ(ierr); 268335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",PETSC_NULL);CHKERRQ(ierr); 269335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C","",PETSC_NULL);CHKERRQ(ierr); 270335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",PETSC_NULL);CHKERRQ(ierr); 271335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",PETSC_NULL);CHKERRQ(ierr); 272335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C","",PETSC_NULL);CHKERRQ(ierr); 273335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C","",PETSC_NULL);CHKERRQ(ierr); 274335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",PETSC_NULL);CHKERRQ(ierr); 2757cb94ee6SHong Zhang PetscFunctionReturn(0); 2767cb94ee6SHong Zhang } 2777cb94ee6SHong Zhang 2787cb94ee6SHong Zhang #undef __FUNCT__ 279214bc6a2SJed Brown #define __FUNCT__ "TSSetUp_Sundials" 280214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts) 2817cb94ee6SHong Zhang { 2827cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2837cb94ee6SHong Zhang PetscErrorCode ierr; 28416016685SHong Zhang PetscInt glosize,locsize,i,flag; 2857cb94ee6SHong Zhang PetscScalar *y_data,*parray; 28616016685SHong Zhang void *mem; 287dcbc6d53SSean Farley PC pc; 28816016685SHong Zhang const PCType pctype; 289ace3abfcSBarry Smith PetscBool pcnone; 2907cb94ee6SHong Zhang 2917cb94ee6SHong Zhang PetscFunctionBegin; 2927cb94ee6SHong Zhang /* get the vector size */ 2937cb94ee6SHong Zhang ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr); 2947cb94ee6SHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr); 2957cb94ee6SHong Zhang 2967cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 297a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 298e32f2f54SBarry Smith if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated"); 2997cb94ee6SHong Zhang 30028adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 3017cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr); 3027cb94ee6SHong Zhang y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y); 3037cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 3047cb94ee6SHong Zhang ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 30542528757SHong Zhang 3067cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr); 3070679a0aeSJed Brown ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr); 3087cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr); 3090679a0aeSJed Brown ierr = PetscLogObjectParent(ts,cvode->ydot);CHKERRQ(ierr); 3107cb94ee6SHong Zhang 3117cb94ee6SHong Zhang /* 3127cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 3137cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 3147cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 3157cb94ee6SHong Zhang */ 3167adad957SLisandro Dalcin ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr); 3177adad957SLisandro Dalcin ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr); 3187cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr); 3197cb94ee6SHong Zhang ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr); 32016016685SHong Zhang 32116016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 32216016685SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 323e32f2f54SBarry Smith if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails"); 32416016685SHong Zhang cvode->mem = mem; 32516016685SHong Zhang 32616016685SHong Zhang /* Set the pointer to user-defined data */ 32716016685SHong Zhang flag = CVodeSetUserData(mem, ts); 328e32f2f54SBarry Smith if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails"); 32916016685SHong Zhang 330fc6b9e64SJed Brown /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 331cdbf8f93SLisandro Dalcin flag = CVodeSetInitStep(mem,(realtype)ts->time_step); 332fc6b9e64SJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed"); 333f1cd61daSJed Brown if (cvode->mindt > 0) { 334f1cd61daSJed Brown flag = CVodeSetMinStep(mem,(realtype)cvode->mindt); 3359f94935aSHong Zhang if (flag){ 3369f94935aSHong Zhang if (flag == CV_MEM_NULL){ 3379f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 3389f94935aSHong Zhang } else if (flag == CV_ILL_INPUT){ 3399f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size"); 3409f94935aSHong Zhang } else { 3419f94935aSHong Zhang SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed"); 3429f94935aSHong Zhang } 3439f94935aSHong Zhang } 344f1cd61daSJed Brown } 345f1cd61daSJed Brown if (cvode->maxdt > 0) { 346f1cd61daSJed Brown flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt); 347f1cd61daSJed Brown if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed"); 348f1cd61daSJed Brown } 349f1cd61daSJed Brown 35016016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 35116016685SHong Zhang * user's right hand side function in u'=f(t,u), the inital time T0, and 35216016685SHong Zhang * the initial dependent variable vector cvode->y */ 35316016685SHong Zhang flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 35416016685SHong Zhang if (flag){ 355e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag); 35616016685SHong Zhang } 35716016685SHong Zhang 3589f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 35916016685SHong Zhang flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 36016016685SHong Zhang if (flag){ 361e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag); 36216016685SHong Zhang } 36316016685SHong Zhang 3649f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 3659f94935aSHong Zhang flag = CVodeSetMaxNumSteps(mem,ts->max_steps); 36616016685SHong Zhang 36716016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 36816016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 369dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr); 370dcbc6d53SSean Farley ierr = PCGetType(pc,&pctype);CHKERRQ(ierr); 371dcbc6d53SSean Farley ierr = PetscTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr); 37216016685SHong Zhang if (pcnone){ 37316016685SHong Zhang flag = CVSpgmr(mem,PREC_NONE,0); 374e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 37516016685SHong Zhang } else { 376f61b2b6aSHong Zhang flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl); 377e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 37816016685SHong Zhang 37916016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 38016016685SHong Zhang and the pointer to the user-defined block data */ 38116016685SHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 382e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 38316016685SHong Zhang } 38416016685SHong Zhang 38516016685SHong Zhang flag = CVSpilsSetGSType(mem, MODIFIED_GS); 386e32f2f54SBarry Smith if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag); 3877cb94ee6SHong Zhang PetscFunctionReturn(0); 3887cb94ee6SHong Zhang } 3897cb94ee6SHong Zhang 3906fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 391dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0}; 3926fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 393dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0}; 394a04cf4d8SBarry Smith 3957cb94ee6SHong Zhang #undef __FUNCT__ 396214bc6a2SJed Brown #define __FUNCT__ "TSSetFromOptions_Sundials" 397214bc6a2SJed Brown PetscErrorCode TSSetFromOptions_Sundials(TS ts) 3987cb94ee6SHong Zhang { 3997cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4007cb94ee6SHong Zhang PetscErrorCode ierr; 4017cb94ee6SHong Zhang int indx; 402ace3abfcSBarry Smith PetscBool flag; 4037bda3a07SJed Brown PC pc; 4047cb94ee6SHong Zhang 4057cb94ee6SHong Zhang PetscFunctionBegin; 4067cb94ee6SHong Zhang ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr); 4076fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr); 4087cb94ee6SHong Zhang if (flag) { 4096fadb2cdSHong Zhang ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr); 4107cb94ee6SHong Zhang } 4116fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr); 4127cb94ee6SHong Zhang if (flag) { 4137cb94ee6SHong Zhang ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 4147cb94ee6SHong Zhang } 4157cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr); 4167cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr); 417f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr); 418f1cd61daSJed Brown ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr); 4197cb94ee6SHong Zhang ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr); 420f61b2b6aSHong Zhang ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,&flag);CHKERRQ(ierr); 421acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr); 4227cb94ee6SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4237bda3a07SJed Brown ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 4247bda3a07SJed Brown ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 4257cb94ee6SHong Zhang PetscFunctionReturn(0); 4267cb94ee6SHong Zhang } 4277cb94ee6SHong Zhang 4287cb94ee6SHong Zhang #undef __FUNCT__ 4297cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials" 4307cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 4317cb94ee6SHong Zhang { 4327cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4337cb94ee6SHong Zhang PetscErrorCode ierr; 4347cb94ee6SHong Zhang char *type; 4357cb94ee6SHong Zhang char atype[] = "Adams"; 4367cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 437ace3abfcSBarry Smith PetscBool iascii,isstring; 4382c823083SHong Zhang long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 4392c823083SHong Zhang PetscInt qlast,qcur; 4405d47aee6SHong Zhang PetscReal hinused,hlast,hcur,tcur,tolsfac; 44142528757SHong Zhang PC pc; 4427cb94ee6SHong Zhang 4437cb94ee6SHong Zhang PetscFunctionBegin; 4447cb94ee6SHong Zhang if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;} 4457cb94ee6SHong Zhang else {type = btype;} 4467cb94ee6SHong Zhang 4472692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4482692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 4497cb94ee6SHong Zhang if (iascii) { 4507cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr); 4517cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr); 4527cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr); 4537cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr); 454f61b2b6aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr); 4557cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 4567cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4577cb94ee6SHong Zhang } else { 4587cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4597cb94ee6SHong Zhang } 460f1cd61daSJed Brown if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);} 461f1cd61daSJed Brown if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);} 4622c823083SHong Zhang 4635d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4645d47aee6SHong Zhang ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr); 4655d47aee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr); 4662c823083SHong Zhang ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals, 4672c823083SHong Zhang &nlinsetups,&nfails,&qlast,&qcur, 4682c823083SHong Zhang &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr); 4692c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr); 4702c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr); 4712c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr); 4722c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr); 4732c823083SHong Zhang 4742c823083SHong Zhang ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr); 4752c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr); 4762c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr); 4772c823083SHong Zhang 4782c823083SHong Zhang ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */ 4792c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr); 4802c823083SHong Zhang ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr); 4812c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr); 48242528757SHong Zhang 48342528757SHong Zhang ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr); 48442528757SHong Zhang ierr = PCView(pc,viewer);CHKERRQ(ierr); 4852c823083SHong Zhang ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4862c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr); 4872c823083SHong Zhang ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr); 4882c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr); 48942528757SHong Zhang 4902c823083SHong Zhang ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4912c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr); 4922c823083SHong Zhang ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4932c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr); 4947cb94ee6SHong Zhang } else if (isstring) { 4957cb94ee6SHong Zhang ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 4967cb94ee6SHong Zhang } 4977cb94ee6SHong Zhang PetscFunctionReturn(0); 4987cb94ee6SHong Zhang } 4997cb94ee6SHong Zhang 5007cb94ee6SHong Zhang 5017cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 5027cb94ee6SHong Zhang EXTERN_C_BEGIN 5037cb94ee6SHong Zhang #undef __FUNCT__ 5047cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials" 5057087cfbeSBarry Smith PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 5067cb94ee6SHong Zhang { 5077cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5087cb94ee6SHong Zhang 5097cb94ee6SHong Zhang PetscFunctionBegin; 5107cb94ee6SHong Zhang cvode->cvode_type = type; 5117cb94ee6SHong Zhang PetscFunctionReturn(0); 5127cb94ee6SHong Zhang } 5137cb94ee6SHong Zhang EXTERN_C_END 5147cb94ee6SHong Zhang 5157cb94ee6SHong Zhang EXTERN_C_BEGIN 5167cb94ee6SHong Zhang #undef __FUNCT__ 517f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl_Sundials" 518f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl) 5197cb94ee6SHong Zhang { 5207cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5217cb94ee6SHong Zhang 5227cb94ee6SHong Zhang PetscFunctionBegin; 523f61b2b6aSHong Zhang cvode->maxl = maxl; 5247cb94ee6SHong Zhang PetscFunctionReturn(0); 5257cb94ee6SHong Zhang } 5267cb94ee6SHong Zhang EXTERN_C_END 5277cb94ee6SHong Zhang 5287cb94ee6SHong Zhang EXTERN_C_BEGIN 5297cb94ee6SHong Zhang #undef __FUNCT__ 5307cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials" 5317087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,double tol) 5327cb94ee6SHong Zhang { 5337cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5347cb94ee6SHong Zhang 5357cb94ee6SHong Zhang PetscFunctionBegin; 5367cb94ee6SHong Zhang cvode->linear_tol = tol; 5377cb94ee6SHong Zhang PetscFunctionReturn(0); 5387cb94ee6SHong Zhang } 5397cb94ee6SHong Zhang EXTERN_C_END 5407cb94ee6SHong Zhang 5417cb94ee6SHong Zhang EXTERN_C_BEGIN 5427cb94ee6SHong Zhang #undef __FUNCT__ 5437cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials" 5447087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 5457cb94ee6SHong Zhang { 5467cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5477cb94ee6SHong Zhang 5487cb94ee6SHong Zhang PetscFunctionBegin; 5497cb94ee6SHong Zhang cvode->gtype = type; 5507cb94ee6SHong Zhang PetscFunctionReturn(0); 5517cb94ee6SHong Zhang } 5527cb94ee6SHong Zhang EXTERN_C_END 5537cb94ee6SHong Zhang 5547cb94ee6SHong Zhang EXTERN_C_BEGIN 5557cb94ee6SHong Zhang #undef __FUNCT__ 5567cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials" 5577087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel) 5587cb94ee6SHong Zhang { 5597cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5607cb94ee6SHong Zhang 5617cb94ee6SHong Zhang PetscFunctionBegin; 5627cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 5637cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 5647cb94ee6SHong Zhang PetscFunctionReturn(0); 5657cb94ee6SHong Zhang } 5667cb94ee6SHong Zhang EXTERN_C_END 5677cb94ee6SHong Zhang 5687cb94ee6SHong Zhang EXTERN_C_BEGIN 5697cb94ee6SHong Zhang #undef __FUNCT__ 570f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials" 5717087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 572f1cd61daSJed Brown { 573f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 574f1cd61daSJed Brown 575f1cd61daSJed Brown PetscFunctionBegin; 576f1cd61daSJed Brown cvode->mindt = mindt; 577f1cd61daSJed Brown PetscFunctionReturn(0); 578f1cd61daSJed Brown } 579f1cd61daSJed Brown EXTERN_C_END 580f1cd61daSJed Brown 581f1cd61daSJed Brown EXTERN_C_BEGIN 582f1cd61daSJed Brown #undef __FUNCT__ 583f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials" 5847087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 585f1cd61daSJed Brown { 586f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 587f1cd61daSJed Brown 588f1cd61daSJed Brown PetscFunctionBegin; 589f1cd61daSJed Brown cvode->maxdt = maxdt; 590f1cd61daSJed Brown PetscFunctionReturn(0); 591f1cd61daSJed Brown } 592f1cd61daSJed Brown EXTERN_C_END 593f1cd61daSJed Brown 594f1cd61daSJed Brown EXTERN_C_BEGIN 595f1cd61daSJed Brown #undef __FUNCT__ 5967cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials" 5977087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc) 5987cb94ee6SHong Zhang { 599dcbc6d53SSean Farley SNES snes; 600dcbc6d53SSean Farley KSP ksp; 601dcbc6d53SSean Farley PetscErrorCode ierr; 6027cb94ee6SHong Zhang 6037cb94ee6SHong Zhang PetscFunctionBegin; 604dcbc6d53SSean Farley ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 605dcbc6d53SSean Farley ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 606dcbc6d53SSean Farley ierr = KSPGetPC(ksp,pc); CHKERRQ(ierr); 6077cb94ee6SHong Zhang PetscFunctionReturn(0); 6087cb94ee6SHong Zhang } 6097cb94ee6SHong Zhang EXTERN_C_END 6107cb94ee6SHong Zhang 6117cb94ee6SHong Zhang EXTERN_C_BEGIN 6127cb94ee6SHong Zhang #undef __FUNCT__ 6137cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials" 6147087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 6157cb94ee6SHong Zhang { 6167cb94ee6SHong Zhang PetscFunctionBegin; 6172c823083SHong Zhang if (nonlin) *nonlin = ts->nonlinear_its; 6182c823083SHong Zhang if (lin) *lin = ts->linear_its; 6197cb94ee6SHong Zhang PetscFunctionReturn(0); 6207cb94ee6SHong Zhang } 6217cb94ee6SHong Zhang EXTERN_C_END 6227cb94ee6SHong Zhang 6237cb94ee6SHong Zhang EXTERN_C_BEGIN 6247cb94ee6SHong Zhang #undef __FUNCT__ 6252bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials" 6267087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s) 6272bfc04deSHong Zhang { 6282bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 6292bfc04deSHong Zhang 6302bfc04deSHong Zhang PetscFunctionBegin; 6312bfc04deSHong Zhang cvode->monitorstep = s; 6322bfc04deSHong Zhang PetscFunctionReturn(0); 6332bfc04deSHong Zhang } 6342bfc04deSHong Zhang EXTERN_C_END 6357cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 6367cb94ee6SHong Zhang 6377cb94ee6SHong Zhang #undef __FUNCT__ 6387cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations" 6397cb94ee6SHong Zhang /*@C 6407cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 6417cb94ee6SHong Zhang 6427cb94ee6SHong Zhang Not Collective 6437cb94ee6SHong Zhang 6447cb94ee6SHong Zhang Input parameters: 6457cb94ee6SHong Zhang . ts - the time-step context 6467cb94ee6SHong Zhang 6477cb94ee6SHong Zhang Output Parameters: 6487cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 6497cb94ee6SHong Zhang - lin - number of linear iterations 6507cb94ee6SHong Zhang 6517cb94ee6SHong Zhang Level: advanced 6527cb94ee6SHong Zhang 6537cb94ee6SHong Zhang Notes: 6547cb94ee6SHong Zhang These return the number since the creation of the TS object 6557cb94ee6SHong Zhang 6567cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations 6577cb94ee6SHong Zhang 658f61b2b6aSHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetMaxl(), 6597cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 660f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 661a43b19c4SJed Brown TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime() 6627cb94ee6SHong Zhang 6637cb94ee6SHong Zhang @*/ 6647087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 6657cb94ee6SHong Zhang { 6664ac538c5SBarry Smith PetscErrorCode ierr; 6677cb94ee6SHong Zhang 6687cb94ee6SHong Zhang PetscFunctionBegin; 6694ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr); 6707cb94ee6SHong Zhang PetscFunctionReturn(0); 6717cb94ee6SHong Zhang } 6727cb94ee6SHong Zhang 6737cb94ee6SHong Zhang #undef __FUNCT__ 6747cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType" 6757cb94ee6SHong Zhang /*@ 6767cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 6777cb94ee6SHong Zhang 6783f9fe445SBarry Smith Logically Collective on TS 6797cb94ee6SHong Zhang 6807cb94ee6SHong Zhang Input parameters: 6817cb94ee6SHong Zhang + ts - the time-step context 6827cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 6837cb94ee6SHong Zhang 6847cb94ee6SHong Zhang Level: intermediate 6857cb94ee6SHong Zhang 6867cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula 6877cb94ee6SHong Zhang 688f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetMaxl(), 6897cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 690f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 6917cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 692a43b19c4SJed Brown TSSetExactFinalTime() 6937cb94ee6SHong Zhang @*/ 6947087cfbeSBarry Smith PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type) 6957cb94ee6SHong Zhang { 6964ac538c5SBarry Smith PetscErrorCode ierr; 6977cb94ee6SHong Zhang 6987cb94ee6SHong Zhang PetscFunctionBegin; 6994ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr); 7007cb94ee6SHong Zhang PetscFunctionReturn(0); 7017cb94ee6SHong Zhang } 7027cb94ee6SHong Zhang 7037cb94ee6SHong Zhang #undef __FUNCT__ 704f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl" 7057cb94ee6SHong Zhang /*@ 706f61b2b6aSHong Zhang TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 7077cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 708f61b2b6aSHong Zhang this is the maximum number of GMRES steps that will be used. 7097cb94ee6SHong Zhang 7103f9fe445SBarry Smith Logically Collective on TS 7117cb94ee6SHong Zhang 7127cb94ee6SHong Zhang Input parameters: 7137cb94ee6SHong Zhang + ts - the time-step context 714f61b2b6aSHong Zhang - maxl - number of direction vectors (the dimension of Krylov subspace). 7157cb94ee6SHong Zhang 7167cb94ee6SHong Zhang Level: advanced 7177cb94ee6SHong Zhang 718f61b2b6aSHong Zhang .keywords: GMRES 7197cb94ee6SHong Zhang 7207cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 7217cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 722f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7237cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 724a43b19c4SJed Brown TSSetExactFinalTime() 7257cb94ee6SHong Zhang 7267cb94ee6SHong Zhang @*/ 727f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl(TS ts,PetscInt maxl) 7287cb94ee6SHong Zhang { 7294ac538c5SBarry Smith PetscErrorCode ierr; 7307cb94ee6SHong Zhang 7317cb94ee6SHong Zhang PetscFunctionBegin; 732f61b2b6aSHong Zhang PetscValidLogicalCollectiveInt(ts,maxl,2); 733f61b2b6aSHong Zhang ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));CHKERRQ(ierr); 7347cb94ee6SHong Zhang PetscFunctionReturn(0); 7357cb94ee6SHong Zhang } 7367cb94ee6SHong Zhang 7377cb94ee6SHong Zhang #undef __FUNCT__ 7387cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance" 7397cb94ee6SHong Zhang /*@ 7407cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 7417cb94ee6SHong Zhang system by SUNDIALS. 7427cb94ee6SHong Zhang 7433f9fe445SBarry Smith Logically Collective on TS 7447cb94ee6SHong Zhang 7457cb94ee6SHong Zhang Input parameters: 7467cb94ee6SHong Zhang + ts - the time-step context 7477cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 7487cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 7497cb94ee6SHong Zhang 7507cb94ee6SHong Zhang Level: advanced 7517cb94ee6SHong Zhang 7527cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS 7537cb94ee6SHong Zhang 754f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7557cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 756f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7577cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 758a43b19c4SJed Brown TSSetExactFinalTime() 7597cb94ee6SHong Zhang 7607cb94ee6SHong Zhang @*/ 7617087cfbeSBarry Smith PetscErrorCode TSSundialsSetLinearTolerance(TS ts,double tol) 7627cb94ee6SHong Zhang { 7634ac538c5SBarry Smith PetscErrorCode ierr; 7647cb94ee6SHong Zhang 7657cb94ee6SHong Zhang PetscFunctionBegin; 766c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,tol,2); 7674ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr); 7687cb94ee6SHong Zhang PetscFunctionReturn(0); 7697cb94ee6SHong Zhang } 7707cb94ee6SHong Zhang 7717cb94ee6SHong Zhang #undef __FUNCT__ 7727cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType" 7737cb94ee6SHong Zhang /*@ 7747cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 7757cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 7767cb94ee6SHong Zhang 7773f9fe445SBarry Smith Logically Collective on TS 7787cb94ee6SHong Zhang 7797cb94ee6SHong Zhang Input parameters: 7807cb94ee6SHong Zhang + ts - the time-step context 7817cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 7827cb94ee6SHong Zhang 7837cb94ee6SHong Zhang Level: advanced 7847cb94ee6SHong Zhang 7857cb94ee6SHong Zhang .keywords: Sundials, orthogonalization 7867cb94ee6SHong Zhang 787f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7887cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 789f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7907cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 791a43b19c4SJed Brown TSSetExactFinalTime() 7927cb94ee6SHong Zhang 7937cb94ee6SHong Zhang @*/ 7947087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 7957cb94ee6SHong Zhang { 7964ac538c5SBarry Smith PetscErrorCode ierr; 7977cb94ee6SHong Zhang 7987cb94ee6SHong Zhang PetscFunctionBegin; 7994ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr); 8007cb94ee6SHong Zhang PetscFunctionReturn(0); 8017cb94ee6SHong Zhang } 8027cb94ee6SHong Zhang 8037cb94ee6SHong Zhang #undef __FUNCT__ 8047cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance" 8057cb94ee6SHong Zhang /*@ 8067cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 8077cb94ee6SHong Zhang Sundials for error control. 8087cb94ee6SHong Zhang 8093f9fe445SBarry Smith Logically Collective on TS 8107cb94ee6SHong Zhang 8117cb94ee6SHong Zhang Input parameters: 8127cb94ee6SHong Zhang + ts - the time-step context 8137cb94ee6SHong Zhang . aabs - the absolute tolerance 8147cb94ee6SHong Zhang - rel - the relative tolerance 8157cb94ee6SHong Zhang 8167cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 8177cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 8187cb94ee6SHong Zhang 8197cb94ee6SHong Zhang Level: intermediate 8207cb94ee6SHong Zhang 8217cb94ee6SHong Zhang .keywords: Sundials, tolerance 8227cb94ee6SHong Zhang 823f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(), 8247cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 825f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8267cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 827a43b19c4SJed Brown TSSetExactFinalTime() 8287cb94ee6SHong Zhang 8297cb94ee6SHong Zhang @*/ 8307087cfbeSBarry Smith PetscErrorCode TSSundialsSetTolerance(TS ts,double aabs,double rel) 8317cb94ee6SHong Zhang { 8324ac538c5SBarry Smith PetscErrorCode ierr; 8337cb94ee6SHong Zhang 8347cb94ee6SHong Zhang PetscFunctionBegin; 8354ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr); 8367cb94ee6SHong Zhang PetscFunctionReturn(0); 8377cb94ee6SHong Zhang } 8387cb94ee6SHong Zhang 8397cb94ee6SHong Zhang #undef __FUNCT__ 8407cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC" 8417cb94ee6SHong Zhang /*@ 8427cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 8437cb94ee6SHong Zhang 8447cb94ee6SHong Zhang Input Parameter: 8457cb94ee6SHong Zhang . ts - the time-step context 8467cb94ee6SHong Zhang 8477cb94ee6SHong Zhang Output Parameter: 8487cb94ee6SHong Zhang . pc - the preconditioner context 8497cb94ee6SHong Zhang 8507cb94ee6SHong Zhang Level: advanced 8517cb94ee6SHong Zhang 852f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 8537cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 854f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8557cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 8567cb94ee6SHong Zhang @*/ 8577087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC(TS ts,PC *pc) 8587cb94ee6SHong Zhang { 8594ac538c5SBarry Smith PetscErrorCode ierr; 8607cb94ee6SHong Zhang 8617cb94ee6SHong Zhang PetscFunctionBegin; 8624ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr); 8637cb94ee6SHong Zhang PetscFunctionReturn(0); 8647cb94ee6SHong Zhang } 8657cb94ee6SHong Zhang 8667cb94ee6SHong Zhang #undef __FUNCT__ 867f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep" 868f1cd61daSJed Brown /*@ 869f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 870f1cd61daSJed Brown 871f1cd61daSJed Brown Input Parameter: 872f1cd61daSJed Brown + ts - the time-step context 873f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 874f1cd61daSJed Brown 875fc6b9e64SJed Brown Note: 876fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 877fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 878fc6b9e64SJed Brown 879f1cd61daSJed Brown Level: beginner 880f1cd61daSJed Brown 881f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 882f1cd61daSJed Brown @*/ 8837087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 884f1cd61daSJed Brown { 8854ac538c5SBarry Smith PetscErrorCode ierr; 886f1cd61daSJed Brown 887f1cd61daSJed Brown PetscFunctionBegin; 8884ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr); 889f1cd61daSJed Brown PetscFunctionReturn(0); 890f1cd61daSJed Brown } 891f1cd61daSJed Brown 892f1cd61daSJed Brown #undef __FUNCT__ 893f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep" 894f1cd61daSJed Brown /*@ 895f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 896f1cd61daSJed Brown 897f1cd61daSJed Brown Input Parameter: 898f1cd61daSJed Brown + ts - the time-step context 899f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 900f1cd61daSJed Brown 901f1cd61daSJed Brown Level: beginner 902f1cd61daSJed Brown 903f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 904f1cd61daSJed Brown @*/ 9057087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 906f1cd61daSJed Brown { 9074ac538c5SBarry Smith PetscErrorCode ierr; 908f1cd61daSJed Brown 909f1cd61daSJed Brown PetscFunctionBegin; 9104ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 911f1cd61daSJed Brown PetscFunctionReturn(0); 912f1cd61daSJed Brown } 913f1cd61daSJed Brown 914f1cd61daSJed Brown #undef __FUNCT__ 9152bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps" 9162bfc04deSHong Zhang /*@ 9172bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 9182bfc04deSHong Zhang 9192bfc04deSHong Zhang Input Parameter: 9202bfc04deSHong Zhang + ts - the time-step context 9212bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 9222bfc04deSHong Zhang 9232bfc04deSHong Zhang Level: beginner 9242bfc04deSHong Zhang 925f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 9262bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 927f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 9282bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 9292bfc04deSHong Zhang @*/ 9307087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft) 9312bfc04deSHong Zhang { 9324ac538c5SBarry Smith PetscErrorCode ierr; 9332bfc04deSHong Zhang 9342bfc04deSHong Zhang PetscFunctionBegin; 9354ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 9362bfc04deSHong Zhang PetscFunctionReturn(0); 9372bfc04deSHong Zhang } 9387cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 9397cb94ee6SHong Zhang /*MC 94096f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 9417cb94ee6SHong Zhang 9427cb94ee6SHong Zhang Options Database: 9437cb94ee6SHong Zhang + -ts_sundials_type <bdf,adams> 9447cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 9457cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 9467cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 9477cb94ee6SHong Zhang . -ts_sundials_linear_tolerance <tol> 948f61b2b6aSHong Zhang . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 94916016685SHong Zhang - -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps 95016016685SHong Zhang 9517cb94ee6SHong Zhang 9527cb94ee6SHong Zhang Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply 9537cb94ee6SHong Zhang only PETSc PC options 9547cb94ee6SHong Zhang 9557cb94ee6SHong Zhang Level: beginner 9567cb94ee6SHong Zhang 957f61b2b6aSHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(), 958a43b19c4SJed Brown TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime() 9597cb94ee6SHong Zhang 9607cb94ee6SHong Zhang M*/ 9617cb94ee6SHong Zhang EXTERN_C_BEGIN 9627cb94ee6SHong Zhang #undef __FUNCT__ 9637cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials" 9647087cfbeSBarry Smith PetscErrorCode TSCreate_Sundials(TS ts) 9657cb94ee6SHong Zhang { 9667cb94ee6SHong Zhang TS_Sundials *cvode; 9677cb94ee6SHong Zhang PetscErrorCode ierr; 96842528757SHong Zhang PC pc; 9697cb94ee6SHong Zhang 9707cb94ee6SHong Zhang PetscFunctionBegin; 971277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Sundials; 97228adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 97328adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 974214bc6a2SJed Brown ts->ops->setup = TSSetUp_Sundials; 975b4eba00bSSean Farley ts->ops->step = TSStep_Sundials; 976b4eba00bSSean Farley ts->ops->interpolate = TSInterpolate_Sundials; 977214bc6a2SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Sundials; 9787cb94ee6SHong Zhang 97938f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr); 9807cb94ee6SHong Zhang ts->data = (void*)cvode; 9816fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 9826fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 983f61b2b6aSHong Zhang cvode->maxl = 5; 9847cb94ee6SHong Zhang cvode->linear_tol = .05; 9857cb94ee6SHong Zhang 986b4eba00bSSean Farley cvode->monitorstep = PETSC_TRUE; 9877cb94ee6SHong Zhang 988a07356b8SSatish Balay ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr); 989f1cd61daSJed Brown 990f1cd61daSJed Brown cvode->mindt = -1.; 991f1cd61daSJed Brown cvode->maxdt = -1.; 992f1cd61daSJed Brown 9937cb94ee6SHong Zhang /* set tolerance for Sundials */ 9947cb94ee6SHong Zhang cvode->reltol = 1e-6; 9952c823083SHong Zhang cvode->abstol = 1e-6; 9967cb94ee6SHong Zhang 99742528757SHong Zhang /* set PCNONE as default pctype */ 99842528757SHong Zhang ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr); 99942528757SHong Zhang ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 100042528757SHong Zhang 1001b4eba00bSSean Farley if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE; 1002a43b19c4SJed Brown 10037cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials", 10047cb94ee6SHong Zhang TSSundialsSetType_Sundials);CHKERRQ(ierr); 1005f61b2b6aSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C", 1006f61b2b6aSHong Zhang "TSSundialsSetMaxl_Sundials", 1007f61b2b6aSHong Zhang TSSundialsSetMaxl_Sundials);CHKERRQ(ierr); 10087cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C", 10097cb94ee6SHong Zhang "TSSundialsSetLinearTolerance_Sundials", 10107cb94ee6SHong Zhang TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 10117cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C", 10127cb94ee6SHong Zhang "TSSundialsSetGramSchmidtType_Sundials", 10137cb94ee6SHong Zhang TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 10147cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C", 10157cb94ee6SHong Zhang "TSSundialsSetTolerance_Sundials", 10167cb94ee6SHong Zhang TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 1017f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C", 1018f1cd61daSJed Brown "TSSundialsSetMinTimeStep_Sundials", 1019f1cd61daSJed Brown TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr); 1020f1cd61daSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C", 1021f1cd61daSJed Brown "TSSundialsSetMaxTimeStep_Sundials", 1022f1cd61daSJed Brown TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr); 10237cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C", 10247cb94ee6SHong Zhang "TSSundialsGetPC_Sundials", 10257cb94ee6SHong Zhang TSSundialsGetPC_Sundials);CHKERRQ(ierr); 10267cb94ee6SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C", 10277cb94ee6SHong Zhang "TSSundialsGetIterations_Sundials", 10287cb94ee6SHong Zhang TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 10292bfc04deSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C", 10302bfc04deSHong Zhang "TSSundialsMonitorInternalSteps_Sundials", 10312bfc04deSHong Zhang TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr); 10322bfc04deSHong Zhang 10337cb94ee6SHong Zhang PetscFunctionReturn(0); 10347cb94ee6SHong Zhang } 10357cb94ee6SHong Zhang EXTERN_C_END 1036