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 */ 81c7d2463SJed Brown #include <../src/ts/impls/implicit/sundials/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 */ 14bbd56ea5SKarl Rupp PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,booleantype jok,booleantype *jcurPtr, 15bbd56ea5SKarl Rupp realtype _gamma,void *P_data,N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3) 167cb94ee6SHong Zhang { 177cb94ee6SHong Zhang TS ts = (TS) P_data; 187cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 19dcbc6d53SSean Farley PC pc; 200679a0aeSJed Brown Mat J,P; 210679a0aeSJed Brown Vec yy = cvode->w1,yydot = cvode->ydot; 220679a0aeSJed Brown PetscReal gm = (PetscReal)_gamma; 237cb94ee6SHong Zhang PetscScalar *y_data; 247cb94ee6SHong Zhang 257cb94ee6SHong Zhang PetscFunctionBegin; 265f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetIJacobian(ts,&J,&P,NULL,NULL)); 277cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(y); 285f80ce2aSJacob Faibussowitsch CHKERRQ(VecPlaceArray(yy,y_data)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(yydot)); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */ 300679a0aeSJed Brown /* compute the shifted Jacobian (1/gm)*I + Jrest */ 315f80ce2aSJacob Faibussowitsch CHKERRQ(TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,J,P,PETSC_FALSE)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(VecResetArray(yy)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(P,gm)); /* turn into I-gm*Jrest, J is not used by Sundials */ 347cb94ee6SHong Zhang *jcurPtr = TRUE; 355f80ce2aSJacob Faibussowitsch CHKERRQ(TSSundialsGetPC(ts,&pc)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetOperators(pc,J,P)); 377cb94ee6SHong Zhang PetscFunctionReturn(0); 387cb94ee6SHong Zhang } 397cb94ee6SHong Zhang 407cb94ee6SHong Zhang /* 417cb94ee6SHong Zhang TSPSolve_Sundials - routine that we provide to Sundials that applies the preconditioner. 427cb94ee6SHong Zhang */ 434ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z, 444ac7836bSHong Zhang realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp) 457cb94ee6SHong Zhang { 467cb94ee6SHong Zhang TS ts = (TS) P_data; 477cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 48dcbc6d53SSean Farley PC pc; 497cb94ee6SHong Zhang Vec rr = cvode->w1,zz = cvode->w2; 507cb94ee6SHong Zhang PetscScalar *r_data,*z_data; 517cb94ee6SHong Zhang 527cb94ee6SHong Zhang PetscFunctionBegin; 537cb94ee6SHong Zhang /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/ 547cb94ee6SHong Zhang r_data = (PetscScalar*) N_VGetArrayPointer(r); 557cb94ee6SHong Zhang z_data = (PetscScalar*) N_VGetArrayPointer(z); 565f80ce2aSJacob Faibussowitsch CHKERRQ(VecPlaceArray(rr,r_data)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(VecPlaceArray(zz,z_data)); 584ac7836bSHong Zhang 597cb94ee6SHong Zhang /* Solve the Px=r and put the result in zz */ 605f80ce2aSJacob Faibussowitsch CHKERRQ(TSSundialsGetPC(ts,&pc)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(PCApply(pc,rr,zz)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(VecResetArray(rr)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(VecResetArray(zz)); 647cb94ee6SHong Zhang PetscFunctionReturn(0); 657cb94ee6SHong Zhang } 667cb94ee6SHong Zhang 677cb94ee6SHong Zhang /* 687cb94ee6SHong Zhang TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side. 697cb94ee6SHong Zhang */ 704ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx) 717cb94ee6SHong Zhang { 727cb94ee6SHong Zhang TS ts = (TS) ctx; 735fcef5e4SJed Brown DM dm; 74942e3340SBarry Smith DMTS tsdm; 755fcef5e4SJed Brown TSIFunction ifunction; 76ce94432eSBarry Smith MPI_Comm comm; 777cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 780679a0aeSJed Brown Vec yy = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot; 797cb94ee6SHong Zhang PetscScalar *y_data,*ydot_data; 807cb94ee6SHong Zhang 817cb94ee6SHong Zhang PetscFunctionBegin; 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)ts,&comm)); 835ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 847cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(y); 857cb94ee6SHong Zhang ydot_data = (PetscScalar*) N_VGetArrayPointer(ydot); 865f80ce2aSJacob Faibussowitsch CHKERRABORT(comm,VecPlaceArray(yy,y_data)); 875f80ce2aSJacob Faibussowitsch CHKERRABORT(comm,VecPlaceArray(yyd,ydot_data)); 884ac7836bSHong Zhang 895fcef5e4SJed Brown /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */ 905f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&dm)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDMTS(dm,&tsdm)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSGetIFunction(dm,&ifunction,NULL)); 935fcef5e4SJed Brown if (!ifunction) { 945f80ce2aSJacob Faibussowitsch CHKERRQ(TSComputeRHSFunction(ts,t,yy,yyd)); 95bc0cc02bSJed Brown } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */ 965f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(yydot)); 975f80ce2aSJacob Faibussowitsch CHKERRABORT(comm,TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(yyd,-1.)); 99bc0cc02bSJed Brown } 1005f80ce2aSJacob Faibussowitsch CHKERRABORT(comm,VecResetArray(yy)); 1015f80ce2aSJacob Faibussowitsch CHKERRABORT(comm,VecResetArray(yyd)); 1024ac7836bSHong Zhang PetscFunctionReturn(0); 1037cb94ee6SHong Zhang } 1047cb94ee6SHong Zhang 1057cb94ee6SHong Zhang /* 106b4eba00bSSean Farley TSStep_Sundials - Calls Sundials to integrate the ODE. 1077cb94ee6SHong Zhang */ 108b4eba00bSSean Farley PetscErrorCode TSStep_Sundials(TS ts) 1097cb94ee6SHong Zhang { 1107cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 111b4eba00bSSean Farley PetscInt flag; 112be5899b3SLisandro Dalcin long int nits,lits,nsteps; 1137cb94ee6SHong Zhang realtype t,tout; 1147cb94ee6SHong Zhang PetscScalar *y_data; 1157cb94ee6SHong Zhang void *mem; 1167cb94ee6SHong Zhang 1177cb94ee6SHong Zhang PetscFunctionBegin; 11816016685SHong Zhang mem = cvode->mem; 1197cb94ee6SHong Zhang tout = ts->max_time; 1205f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(ts->vec_sol,&y_data)); 1217cb94ee6SHong Zhang N_VSetArrayPointer((realtype*)y_data,cvode->y); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(ts->vec_sol,NULL)); 123186e87acSLisandro Dalcin 1249687d888SLisandro Dalcin /* We would like to TSPreStage() and TSPostStage() 1259687d888SLisandro Dalcin * before each stage solve but CVode does not appear to support this. */ 126be5899b3SLisandro Dalcin if (cvode->monitorstep) 127be5899b3SLisandro Dalcin flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP); 128be5899b3SLisandro Dalcin else 129be5899b3SLisandro Dalcin flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL); 1309f94935aSHong Zhang 1319f94935aSHong Zhang if (flag) { /* display error message */ 1329f94935aSHong Zhang switch (flag) { 1339f94935aSHong Zhang case CV_ILL_INPUT: 1349f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT"); 1359f94935aSHong Zhang break; 1369f94935aSHong Zhang case CV_TOO_CLOSE: 1379f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE"); 1389f94935aSHong Zhang break; 1393c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1409f94935aSHong Zhang PetscReal tcur; 1415f80ce2aSJacob Faibussowitsch CHKERRQ(CVodeGetNumSteps(mem,&nsteps)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(CVodeGetCurrentTime(mem,&tcur)); 14398921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_WORK. At t=%g, nsteps %D exceeds maxstep %D. Increase '-ts_max_steps <>' or modify TSSetMaxSteps()",(double)tcur,nsteps,ts->max_steps); 1443c7fefeeSJed Brown } break; 1459f94935aSHong Zhang case CV_TOO_MUCH_ACC: 1469f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC"); 1479f94935aSHong Zhang break; 1489f94935aSHong Zhang case CV_ERR_FAILURE: 1499f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE"); 1509f94935aSHong Zhang break; 1519f94935aSHong Zhang case CV_CONV_FAILURE: 1529f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE"); 1539f94935aSHong Zhang break; 1549f94935aSHong Zhang case CV_LINIT_FAIL: 1559f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL"); 1569f94935aSHong Zhang break; 1579f94935aSHong Zhang case CV_LSETUP_FAIL: 1589f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL"); 1599f94935aSHong Zhang break; 1609f94935aSHong Zhang case CV_LSOLVE_FAIL: 1619f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL"); 1629f94935aSHong Zhang break; 1639f94935aSHong Zhang case CV_RHSFUNC_FAIL: 1649f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL"); 1659f94935aSHong Zhang break; 1669f94935aSHong Zhang case CV_FIRST_RHSFUNC_ERR: 1679f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR"); 1689f94935aSHong Zhang break; 1699f94935aSHong Zhang case CV_REPTD_RHSFUNC_ERR: 1709f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR"); 1719f94935aSHong Zhang break; 1729f94935aSHong Zhang case CV_UNREC_RHSFUNC_ERR: 1739f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR"); 1749f94935aSHong Zhang break; 1759f94935aSHong Zhang case CV_RTFUNC_FAIL: 1769f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL"); 1779f94935aSHong Zhang break; 1789f94935aSHong Zhang default: 17998921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag); 1809f94935aSHong Zhang } 1819f94935aSHong Zhang } 1829f94935aSHong Zhang 183be5899b3SLisandro Dalcin /* log inner nonlinear and linear iterations */ 1845f80ce2aSJacob Faibussowitsch CHKERRQ(CVodeGetNumNonlinSolvIters(mem,&nits)); 1855f80ce2aSJacob Faibussowitsch CHKERRQ(CVSpilsGetNumLinIters(mem,&lits)); 186be5899b3SLisandro Dalcin ts->snes_its += nits; ts->ksp_its = lits; 187be5899b3SLisandro Dalcin 1887cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */ 1895f80ce2aSJacob Faibussowitsch CHKERRQ(VecPlaceArray(cvode->w1,y_data)); 1905f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(cvode->w1,cvode->update)); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(VecResetArray(cvode->w1)); 1925f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(cvode->update,ts->vec_sol)); 193186e87acSLisandro Dalcin 194186e87acSLisandro Dalcin ts->time_step = t - ts->ptime; 195186e87acSLisandro Dalcin ts->ptime = t; 196186e87acSLisandro Dalcin 1975f80ce2aSJacob Faibussowitsch CHKERRQ(CVodeGetNumSteps(mem,&nsteps)); 1982808aa04SLisandro Dalcin if (!cvode->monitorstep) ts->steps += nsteps - 1; /* TSStep() increments the step counter by one */ 199b4eba00bSSean Farley PetscFunctionReturn(0); 200b4eba00bSSean Farley } 201b4eba00bSSean Farley 202b4eba00bSSean Farley static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X) 203b4eba00bSSean Farley { 204b4eba00bSSean Farley TS_Sundials *cvode = (TS_Sundials*)ts->data; 205b4eba00bSSean Farley N_Vector y; 206b4eba00bSSean Farley PetscScalar *x_data; 207b4eba00bSSean Farley PetscInt glosize,locsize; 208b4eba00bSSean Farley 209b4eba00bSSean Farley PetscFunctionBegin; 210b4eba00bSSean Farley /* get the vector size */ 2115f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(X,&glosize)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(X,&locsize)); 2135f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(X,&x_data)); 214b4eba00bSSean Farley 2154d78c9e0SClaas Abert /* Initialize N_Vec y with x_data */ 216918687b8SPatrick Sanan if (cvode->use_dense) { 217918687b8SPatrick Sanan PetscMPIInt size; 218918687b8SPatrick Sanan 2195f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 2203c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"TSSUNDIALS only supports a dense solve in the serial case"); 221918687b8SPatrick Sanan y = N_VMake_Serial(locsize,(realtype*)x_data); 222918687b8SPatrick Sanan } else { 2234d78c9e0SClaas Abert y = N_VMake_Parallel(cvode->comm_sundials,locsize,glosize,(realtype*)x_data); 224918687b8SPatrick Sanan } 225918687b8SPatrick Sanan 2263c633725SBarry Smith PetscCheck(y,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Interpolated y is not allocated"); 227b4eba00bSSean Farley 2285f80ce2aSJacob Faibussowitsch CHKERRQ(CVodeGetDky(cvode->mem,t,0,y)); 2295f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(X,&x_data)); 2307cb94ee6SHong Zhang PetscFunctionReturn(0); 2317cb94ee6SHong Zhang } 2327cb94ee6SHong Zhang 233277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts) 234277b19d0SLisandro Dalcin { 235277b19d0SLisandro Dalcin TS_Sundials *cvode = (TS_Sundials*)ts->data; 236277b19d0SLisandro Dalcin 237277b19d0SLisandro Dalcin PetscFunctionBegin; 2385f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&cvode->update)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&cvode->ydot)); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&cvode->w1)); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&cvode->w2)); 242bbd56ea5SKarl Rupp if (cvode->mem) CVodeFree(&cvode->mem); 243277b19d0SLisandro Dalcin PetscFunctionReturn(0); 244277b19d0SLisandro Dalcin } 245277b19d0SLisandro Dalcin 2467cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts) 2477cb94ee6SHong Zhang { 2487cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2497cb94ee6SHong Zhang 2507cb94ee6SHong Zhang PetscFunctionBegin; 2515f80ce2aSJacob Faibussowitsch CHKERRQ(TSReset_Sundials(ts)); 2525f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_free(&(cvode->comm_sundials))); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ts->data)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",NULL)); 2555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",NULL)); 2565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",NULL)); 2575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",NULL)); 2585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",NULL)); 2595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",NULL)); 2605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",NULL)); 2615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",NULL)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",NULL)); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",NULL)); 2647cb94ee6SHong Zhang PetscFunctionReturn(0); 2657cb94ee6SHong Zhang } 2667cb94ee6SHong Zhang 267214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts) 2687cb94ee6SHong Zhang { 2697cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 27016016685SHong Zhang PetscInt glosize,locsize,i,flag; 2717cb94ee6SHong Zhang PetscScalar *y_data,*parray; 27216016685SHong Zhang void *mem; 273dcbc6d53SSean Farley PC pc; 27419fd82e9SBarry Smith PCType pctype; 275ace3abfcSBarry Smith PetscBool pcnone; 2767cb94ee6SHong Zhang 2777cb94ee6SHong Zhang PetscFunctionBegin; 2782c71b3e2SJacob Faibussowitsch PetscCheckFalse(ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for exact final time option 'MATCHSTEP' when using Sundials"); 279236c45dcSShri Abhyankar 2807cb94ee6SHong Zhang /* get the vector size */ 2815f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(ts->vec_sol,&glosize)); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(ts->vec_sol,&locsize)); 2837cb94ee6SHong Zhang 2847cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 285918687b8SPatrick Sanan if (cvode->use_dense) { 286918687b8SPatrick Sanan PetscMPIInt size; 287918687b8SPatrick Sanan 2885f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 2893c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"TSSUNDIALS only supports a dense solve in the serial case"); 290918687b8SPatrick Sanan cvode->y = N_VNew_Serial(locsize); 291918687b8SPatrick Sanan } else { 292a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 293918687b8SPatrick Sanan } 2943c633725SBarry Smith PetscCheck(cvode->y,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cvode->y is not allocated"); 2957cb94ee6SHong Zhang 29628adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 2975f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(ts->vec_sol,&parray)); 2987cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(cvode->y); 2997cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 3005f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(ts->vec_sol,NULL)); 30142528757SHong Zhang 3025f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(ts->vec_sol,&cvode->update)); 3035f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(ts->vec_sol,&cvode->ydot)); 3045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->update)); 3055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->ydot)); 3067cb94ee6SHong Zhang 3077cb94ee6SHong Zhang /* 3087cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 3097cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 3107cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 3117cb94ee6SHong Zhang */ 3125f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,NULL,&cvode->w1)); 3135f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,NULL,&cvode->w2)); 3145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w1)); 3155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w2)); 31616016685SHong Zhang 31716016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 31816016685SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 3193c633725SBarry Smith PetscCheck(mem,PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails"); 32016016685SHong Zhang cvode->mem = mem; 32116016685SHong Zhang 32216016685SHong Zhang /* Set the pointer to user-defined data */ 32316016685SHong Zhang flag = CVodeSetUserData(mem, ts); 324*28b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails"); 32516016685SHong Zhang 326fc6b9e64SJed Brown /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 327cdbf8f93SLisandro Dalcin flag = CVodeSetInitStep(mem,(realtype)ts->time_step); 328*28b400f6SJacob Faibussowitsch PetscCheck(!flag,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetInitStep() failed"); 329f1cd61daSJed Brown if (cvode->mindt > 0) { 330f1cd61daSJed Brown flag = CVodeSetMinStep(mem,(realtype)cvode->mindt); 3319f94935aSHong Zhang if (flag) { 3322c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag == CV_MEM_NULL,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 3332c71b3e2SJacob Faibussowitsch else PetscCheckFalse(flag == CV_ILL_INPUT,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size"); 334ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed"); 3359f94935aSHong Zhang } 336f1cd61daSJed Brown } 337f1cd61daSJed Brown if (cvode->maxdt > 0) { 338f1cd61daSJed Brown flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt); 339*28b400f6SJacob Faibussowitsch PetscCheck(!flag,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMaxStep() failed"); 340f1cd61daSJed Brown } 341f1cd61daSJed Brown 34216016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 343a5b23f4aSJose E. Roman * user's right hand side function in u'=f(t,u), the initial time T0, and 34416016685SHong Zhang * the initial dependent variable vector cvode->y */ 34516016685SHong Zhang flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 346*28b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag); 34716016685SHong Zhang 3489f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 34916016685SHong Zhang flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 350*28b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag); 35116016685SHong Zhang 352c4e80e11SFlorian /* Specify max order of BDF / ADAMS method */ 353c4e80e11SFlorian if (cvode->maxord != PETSC_DEFAULT) { 354c4e80e11SFlorian flag = CVodeSetMaxOrd(mem,cvode->maxord); 355*28b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxOrd() fails, flag %d",flag); 356c4e80e11SFlorian } 357c4e80e11SFlorian 3589f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 3599f94935aSHong Zhang flag = CVodeSetMaxNumSteps(mem,ts->max_steps); 360*28b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxNumSteps() fails, flag %d",flag); 36116016685SHong Zhang 362918687b8SPatrick Sanan if (cvode->use_dense) { 363918687b8SPatrick Sanan /* call CVDense to use a dense linear solver. */ 364918687b8SPatrick Sanan PetscMPIInt size; 365918687b8SPatrick Sanan 3665f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 3673c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"TSSUNDIALS only supports a dense solve in the serial case"); 368918687b8SPatrick Sanan flag = CVDense(mem,locsize); 369*28b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVDense() fails, flag %d",flag); 370918687b8SPatrick Sanan } else { 37116016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 37216016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 3735f80ce2aSJacob Faibussowitsch CHKERRQ(TSSundialsGetPC(ts,&pc)); 3745f80ce2aSJacob Faibussowitsch CHKERRQ(PCGetType(pc,&pctype)); 3755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone)); 37616016685SHong Zhang if (pcnone) { 37716016685SHong Zhang flag = CVSpgmr(mem,PREC_NONE,0); 378*28b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 37916016685SHong Zhang } else { 380f61b2b6aSHong Zhang flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl); 381*28b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 38216016685SHong Zhang 38316016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 38416016685SHong Zhang and the pointer to the user-defined block data */ 38516016685SHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 386*28b400f6SJacob Faibussowitsch PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 38716016685SHong Zhang } 388918687b8SPatrick Sanan } 3897cb94ee6SHong Zhang PetscFunctionReturn(0); 3907cb94ee6SHong Zhang } 3917cb94ee6SHong Zhang 3926fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 393c793f718SLisandro Dalcin const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",NULL}; 3946fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 395c793f718SLisandro Dalcin const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",NULL}; 396a04cf4d8SBarry Smith 3974416b707SBarry Smith PetscErrorCode TSSetFromOptions_Sundials(PetscOptionItems *PetscOptionsObject,TS ts) 3987cb94ee6SHong Zhang { 3997cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4007cb94ee6SHong Zhang int indx; 401ace3abfcSBarry Smith PetscBool flag; 4027bda3a07SJed Brown PC pc; 4037cb94ee6SHong Zhang 4047cb94ee6SHong Zhang PetscFunctionBegin; 4055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHead(PetscOptionsObject,"SUNDIALS ODE solver options")); 4065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag)); 4077cb94ee6SHong Zhang if (flag) { 4085f80ce2aSJacob Faibussowitsch CHKERRQ(TSSundialsSetType(ts,(TSSundialsLmmType)indx)); 4097cb94ee6SHong Zhang } 4105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag)); 4117cb94ee6SHong Zhang if (flag) { 4125f80ce2aSJacob Faibussowitsch CHKERRQ(TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx)); 4137cb94ee6SHong Zhang } 4145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL)); 4155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL)); 4165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL)); 4175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL)); 4185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,NULL)); 4195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-ts_sundials_maxord","Max Order for BDF/Adams method","TSSundialsSetMaxOrd",cvode->maxord,&cvode->maxord,NULL)); 4205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,NULL)); 4215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internal steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL)); 4225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-ts_sundials_use_dense","Use dense internal solver in SUNDIALS (serial only)","TSSundialsSetUseDense",cvode->use_dense,&cvode->use_dense,NULL)); 4235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsTail()); 4245f80ce2aSJacob Faibussowitsch CHKERRQ(TSSundialsGetPC(ts,&pc)); 4255f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetFromOptions(pc)); 4267cb94ee6SHong Zhang PetscFunctionReturn(0); 4277cb94ee6SHong Zhang } 4287cb94ee6SHong Zhang 4297cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 4307cb94ee6SHong Zhang { 4317cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4327cb94ee6SHong Zhang char *type; 4337cb94ee6SHong Zhang char atype[] = "Adams"; 4347cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 435ace3abfcSBarry Smith PetscBool iascii,isstring; 4362c823083SHong Zhang long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 4372c823083SHong Zhang PetscInt qlast,qcur; 4385d47aee6SHong Zhang PetscReal hinused,hlast,hcur,tcur,tolsfac; 43942528757SHong Zhang PC pc; 4407cb94ee6SHong Zhang 4417cb94ee6SHong Zhang PetscFunctionBegin; 442bbd56ea5SKarl Rupp if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype; 443bbd56ea5SKarl Rupp else type = btype; 4447cb94ee6SHong Zhang 4455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 4465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring)); 4477cb94ee6SHong Zhang if (iascii) { 4485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials integrator does not use SNES!\n")); 4495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials integrator type %s\n",type)); 4505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials maxord %D\n",cvode->maxord)); 4515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",(double)cvode->abstol,(double)cvode->reltol)); 452918687b8SPatrick Sanan if (cvode->use_dense) { 4535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials integrator using a dense linear solve\n")); 454918687b8SPatrick Sanan } else { 4555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",(double)cvode->linear_tol)); 4565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl)); 4577cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 4585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n")); 4597cb94ee6SHong Zhang } else { 4605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n")); 4617cb94ee6SHong Zhang } 462918687b8SPatrick Sanan } 4635f80ce2aSJacob Faibussowitsch if (cvode->mindt > 0) CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",(double)cvode->mindt)); 4645f80ce2aSJacob Faibussowitsch if (cvode->maxdt > 0) CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",(double)cvode->maxdt)); 4652c823083SHong Zhang 4665d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4675f80ce2aSJacob Faibussowitsch CHKERRQ(CVodeGetTolScaleFactor(cvode->mem,&tolsfac)); 4685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac)); 4695f80ce2aSJacob Faibussowitsch CHKERRQ(CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,&nlinsetups,&nfails,&qlast,&qcur,&hinused,&hlast,&hcur,&tcur)); 4705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps)); 4715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals)); 4725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups)); 4735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails)); 4742c823083SHong Zhang 4755f80ce2aSJacob Faibussowitsch CHKERRQ(CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails)); 4765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its)); 4775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails)); 478918687b8SPatrick Sanan if (!cvode->use_dense) { 4795f80ce2aSJacob Faibussowitsch CHKERRQ(CVSpilsGetNumLinIters(cvode->mem, &its)); /* its = no. of calls to TSPrecond_Sundials() */ 4805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its)); 4815f80ce2aSJacob Faibussowitsch CHKERRQ(CVSpilsGetNumConvFails(cvode->mem,&itmp)); 4825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp)); 48342528757SHong Zhang 4845f80ce2aSJacob Faibussowitsch CHKERRQ(TSSundialsGetPC(ts,&pc)); 4855f80ce2aSJacob Faibussowitsch CHKERRQ(PCView(pc,viewer)); 4865f80ce2aSJacob Faibussowitsch CHKERRQ(CVSpilsGetNumPrecEvals(cvode->mem,&itmp)); 4875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp)); 4885f80ce2aSJacob Faibussowitsch CHKERRQ(CVSpilsGetNumPrecSolves(cvode->mem,&itmp)); 4895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp)); 490918687b8SPatrick Sanan } 4915f80ce2aSJacob Faibussowitsch CHKERRQ(CVSpilsGetNumJtimesEvals(cvode->mem,&itmp)); 4925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp)); 4935f80ce2aSJacob Faibussowitsch CHKERRQ(CVSpilsGetNumRhsEvals(cvode->mem,&itmp)); 4945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp)); 4957cb94ee6SHong Zhang } else if (isstring) { 4965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerStringSPrintf(viewer,"Sundials type %s",type)); 4977cb94ee6SHong Zhang } 4987cb94ee6SHong Zhang PetscFunctionReturn(0); 4997cb94ee6SHong Zhang } 5007cb94ee6SHong Zhang 5017cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 5027087cfbeSBarry Smith PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 5037cb94ee6SHong Zhang { 5047cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5057cb94ee6SHong Zhang 5067cb94ee6SHong Zhang PetscFunctionBegin; 5077cb94ee6SHong Zhang cvode->cvode_type = type; 5087cb94ee6SHong Zhang PetscFunctionReturn(0); 5097cb94ee6SHong Zhang } 5107cb94ee6SHong Zhang 511f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl) 5127cb94ee6SHong Zhang { 5137cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5147cb94ee6SHong Zhang 5157cb94ee6SHong Zhang PetscFunctionBegin; 516f61b2b6aSHong Zhang cvode->maxl = maxl; 5177cb94ee6SHong Zhang PetscFunctionReturn(0); 5187cb94ee6SHong Zhang } 5197cb94ee6SHong Zhang 520590ae059SPatrick Sanan PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,PetscReal tol) 5217cb94ee6SHong Zhang { 5227cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5237cb94ee6SHong Zhang 5247cb94ee6SHong Zhang PetscFunctionBegin; 5257cb94ee6SHong Zhang cvode->linear_tol = tol; 5267cb94ee6SHong Zhang PetscFunctionReturn(0); 5277cb94ee6SHong Zhang } 5287cb94ee6SHong Zhang 5297087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 5307cb94ee6SHong Zhang { 5317cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5327cb94ee6SHong Zhang 5337cb94ee6SHong Zhang PetscFunctionBegin; 5347cb94ee6SHong Zhang cvode->gtype = type; 5357cb94ee6SHong Zhang PetscFunctionReturn(0); 5367cb94ee6SHong Zhang } 5377cb94ee6SHong Zhang 538590ae059SPatrick Sanan PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,PetscReal aabs,PetscReal rel) 5397cb94ee6SHong Zhang { 5407cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5417cb94ee6SHong Zhang 5427cb94ee6SHong Zhang PetscFunctionBegin; 5437cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 5447cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 5457cb94ee6SHong Zhang PetscFunctionReturn(0); 5467cb94ee6SHong Zhang } 5477cb94ee6SHong Zhang 5487087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 549f1cd61daSJed Brown { 550f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 551f1cd61daSJed Brown 552f1cd61daSJed Brown PetscFunctionBegin; 553f1cd61daSJed Brown cvode->mindt = mindt; 554f1cd61daSJed Brown PetscFunctionReturn(0); 555f1cd61daSJed Brown } 556f1cd61daSJed Brown 5577087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 558f1cd61daSJed Brown { 559f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 560f1cd61daSJed Brown 561f1cd61daSJed Brown PetscFunctionBegin; 562f1cd61daSJed Brown cvode->maxdt = maxdt; 563f1cd61daSJed Brown PetscFunctionReturn(0); 564f1cd61daSJed Brown } 565918687b8SPatrick Sanan 566918687b8SPatrick Sanan PetscErrorCode TSSundialsSetUseDense_Sundials(TS ts,PetscBool use_dense) 567918687b8SPatrick Sanan { 568918687b8SPatrick Sanan TS_Sundials *cvode = (TS_Sundials*)ts->data; 569918687b8SPatrick Sanan 570918687b8SPatrick Sanan PetscFunctionBegin; 571918687b8SPatrick Sanan cvode->use_dense = use_dense; 572918687b8SPatrick Sanan PetscFunctionReturn(0); 573918687b8SPatrick Sanan } 574918687b8SPatrick Sanan 5757087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc) 5767cb94ee6SHong Zhang { 577dcbc6d53SSean Farley SNES snes; 578dcbc6d53SSean Farley KSP ksp; 5797cb94ee6SHong Zhang 5807cb94ee6SHong Zhang PetscFunctionBegin; 5815f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSNES(ts,&snes)); 5825f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetKSP(snes,&ksp)); 5835f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(ksp,pc)); 5847cb94ee6SHong Zhang PetscFunctionReturn(0); 5857cb94ee6SHong Zhang } 5867cb94ee6SHong Zhang 5877087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 5887cb94ee6SHong Zhang { 5897cb94ee6SHong Zhang PetscFunctionBegin; 5905ef26d82SJed Brown if (nonlin) *nonlin = ts->snes_its; 5915ef26d82SJed Brown if (lin) *lin = ts->ksp_its; 5927cb94ee6SHong Zhang PetscFunctionReturn(0); 5937cb94ee6SHong Zhang } 5947cb94ee6SHong Zhang 5957087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s) 5962bfc04deSHong Zhang { 5972bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5982bfc04deSHong Zhang 5992bfc04deSHong Zhang PetscFunctionBegin; 6002bfc04deSHong Zhang cvode->monitorstep = s; 6012bfc04deSHong Zhang PetscFunctionReturn(0); 6022bfc04deSHong Zhang } 6037cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 6047cb94ee6SHong Zhang 6057cb94ee6SHong Zhang /*@C 6067cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 6077cb94ee6SHong Zhang 6087cb94ee6SHong Zhang Not Collective 6097cb94ee6SHong Zhang 61097bb3fdcSJose E. Roman Input Parameter: 6117cb94ee6SHong Zhang . ts - the time-step context 6127cb94ee6SHong Zhang 6137cb94ee6SHong Zhang Output Parameters: 6147cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 6157cb94ee6SHong Zhang - lin - number of linear iterations 6167cb94ee6SHong Zhang 6177cb94ee6SHong Zhang Level: advanced 6187cb94ee6SHong Zhang 6197cb94ee6SHong Zhang Notes: 6207cb94ee6SHong Zhang These return the number since the creation of the TS object 6217cb94ee6SHong Zhang 622f61b2b6aSHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetMaxl(), 6237cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 624f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 625a43b19c4SJed Brown TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime() 6267cb94ee6SHong Zhang 6277cb94ee6SHong Zhang @*/ 6287087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 6297cb94ee6SHong Zhang { 6307cb94ee6SHong Zhang PetscFunctionBegin; 6315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin))); 6327cb94ee6SHong Zhang PetscFunctionReturn(0); 6337cb94ee6SHong Zhang } 6347cb94ee6SHong Zhang 6357cb94ee6SHong Zhang /*@ 6367cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 6377cb94ee6SHong Zhang 6383f9fe445SBarry Smith Logically Collective on TS 6397cb94ee6SHong Zhang 6408f7d6fe5SPatrick Sanan Input Parameters: 6417cb94ee6SHong Zhang + ts - the time-step context 6427cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 6437cb94ee6SHong Zhang 6447cb94ee6SHong Zhang Level: intermediate 6457cb94ee6SHong Zhang 646f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetMaxl(), 6477cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 648f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 6497cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 650a43b19c4SJed Brown TSSetExactFinalTime() 6517cb94ee6SHong Zhang @*/ 6527087cfbeSBarry Smith PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type) 6537cb94ee6SHong Zhang { 6547cb94ee6SHong Zhang PetscFunctionBegin; 6555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type))); 6567cb94ee6SHong Zhang PetscFunctionReturn(0); 6577cb94ee6SHong Zhang } 6587cb94ee6SHong Zhang 6597cb94ee6SHong Zhang /*@ 660c4e80e11SFlorian TSSundialsSetMaxord - Sets the maximum order for BDF/Adams method used by SUNDIALS. 661c4e80e11SFlorian 662c4e80e11SFlorian Logically Collective on TS 663c4e80e11SFlorian 6648f7d6fe5SPatrick Sanan Input Parameters: 665c4e80e11SFlorian + ts - the time-step context 666c4e80e11SFlorian - maxord - maximum order of BDF / Adams method 667c4e80e11SFlorian 668c4e80e11SFlorian Level: advanced 669c4e80e11SFlorian 670c4e80e11SFlorian .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 671c4e80e11SFlorian TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 672c4e80e11SFlorian TSSundialsGetIterations(), TSSundialsSetType(), 673c4e80e11SFlorian TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 674c4e80e11SFlorian TSSetExactFinalTime() 675c4e80e11SFlorian 676c4e80e11SFlorian @*/ 677c4e80e11SFlorian PetscErrorCode TSSundialsSetMaxord(TS ts,PetscInt maxord) 678c4e80e11SFlorian { 679c4e80e11SFlorian PetscFunctionBegin; 680c4e80e11SFlorian PetscValidLogicalCollectiveInt(ts,maxord,2); 6815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(ts,"TSSundialsSetMaxOrd_C",(TS,PetscInt),(ts,maxord))); 682c4e80e11SFlorian PetscFunctionReturn(0); 683c4e80e11SFlorian } 684c4e80e11SFlorian 685c4e80e11SFlorian /*@ 686f61b2b6aSHong Zhang TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 6877cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 688f61b2b6aSHong Zhang this is the maximum number of GMRES steps that will be used. 6897cb94ee6SHong Zhang 6903f9fe445SBarry Smith Logically Collective on TS 6917cb94ee6SHong Zhang 6928f7d6fe5SPatrick Sanan Input Parameters: 6937cb94ee6SHong Zhang + ts - the time-step context 694f61b2b6aSHong Zhang - maxl - number of direction vectors (the dimension of Krylov subspace). 6957cb94ee6SHong Zhang 6967cb94ee6SHong Zhang Level: advanced 6977cb94ee6SHong Zhang 6987cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 6997cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 700f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7017cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 702a43b19c4SJed Brown TSSetExactFinalTime() 7037cb94ee6SHong Zhang 7047cb94ee6SHong Zhang @*/ 705f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl(TS ts,PetscInt maxl) 7067cb94ee6SHong Zhang { 7077cb94ee6SHong Zhang PetscFunctionBegin; 708f61b2b6aSHong Zhang PetscValidLogicalCollectiveInt(ts,maxl,2); 7095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl))); 7107cb94ee6SHong Zhang PetscFunctionReturn(0); 7117cb94ee6SHong Zhang } 7127cb94ee6SHong Zhang 7137cb94ee6SHong Zhang /*@ 7147cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 7157cb94ee6SHong Zhang system by SUNDIALS. 7167cb94ee6SHong Zhang 7173f9fe445SBarry Smith Logically Collective on TS 7187cb94ee6SHong Zhang 7198f7d6fe5SPatrick Sanan Input Parameters: 7207cb94ee6SHong Zhang + ts - the time-step context 7217cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 7227cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 7237cb94ee6SHong Zhang 7247cb94ee6SHong Zhang Level: advanced 7257cb94ee6SHong Zhang 726f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7277cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 728f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7297cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 730a43b19c4SJed Brown TSSetExactFinalTime() 7317cb94ee6SHong Zhang 7327cb94ee6SHong Zhang @*/ 733590ae059SPatrick Sanan PetscErrorCode TSSundialsSetLinearTolerance(TS ts,PetscReal tol) 7347cb94ee6SHong Zhang { 7357cb94ee6SHong Zhang PetscFunctionBegin; 736c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,tol,2); 7375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,PetscReal),(ts,tol))); 7387cb94ee6SHong Zhang PetscFunctionReturn(0); 7397cb94ee6SHong Zhang } 7407cb94ee6SHong Zhang 7417cb94ee6SHong Zhang /*@ 7427cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 7437cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 7447cb94ee6SHong Zhang 7453f9fe445SBarry Smith Logically Collective on TS 7467cb94ee6SHong Zhang 7478f7d6fe5SPatrick Sanan Input Parameters: 7487cb94ee6SHong Zhang + ts - the time-step context 7497cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 7507cb94ee6SHong Zhang 7517cb94ee6SHong Zhang Level: advanced 7527cb94ee6SHong Zhang 753f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7547cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 755f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7567cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 757a43b19c4SJed Brown TSSetExactFinalTime() 7587cb94ee6SHong Zhang 7597cb94ee6SHong Zhang @*/ 7607087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 7617cb94ee6SHong Zhang { 7627cb94ee6SHong Zhang PetscFunctionBegin; 7635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type))); 7647cb94ee6SHong Zhang PetscFunctionReturn(0); 7657cb94ee6SHong Zhang } 7667cb94ee6SHong Zhang 7677cb94ee6SHong Zhang /*@ 7687cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 7697cb94ee6SHong Zhang Sundials for error control. 7707cb94ee6SHong Zhang 7713f9fe445SBarry Smith Logically Collective on TS 7727cb94ee6SHong Zhang 7738f7d6fe5SPatrick Sanan Input Parameters: 7747cb94ee6SHong Zhang + ts - the time-step context 7757cb94ee6SHong Zhang . aabs - the absolute tolerance 7767cb94ee6SHong Zhang - rel - the relative tolerance 7777cb94ee6SHong Zhang 7787cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 7797cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 7807cb94ee6SHong Zhang 7817cb94ee6SHong Zhang Level: intermediate 7827cb94ee6SHong Zhang 783f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(), 7847cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 785f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7867cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 787a43b19c4SJed Brown TSSetExactFinalTime() 7887cb94ee6SHong Zhang 7897cb94ee6SHong Zhang @*/ 790590ae059SPatrick Sanan PetscErrorCode TSSundialsSetTolerance(TS ts,PetscReal aabs,PetscReal rel) 7917cb94ee6SHong Zhang { 7927cb94ee6SHong Zhang PetscFunctionBegin; 7935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,PetscReal,PetscReal),(ts,aabs,rel))); 7947cb94ee6SHong Zhang PetscFunctionReturn(0); 7957cb94ee6SHong Zhang } 7967cb94ee6SHong Zhang 7977cb94ee6SHong Zhang /*@ 7987cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 7997cb94ee6SHong Zhang 8007cb94ee6SHong Zhang Input Parameter: 8017cb94ee6SHong Zhang . ts - the time-step context 8027cb94ee6SHong Zhang 8037cb94ee6SHong Zhang Output Parameter: 8047cb94ee6SHong Zhang . pc - the preconditioner context 8057cb94ee6SHong Zhang 8067cb94ee6SHong Zhang Level: advanced 8077cb94ee6SHong Zhang 808f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 8097cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 810f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8117cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 8127cb94ee6SHong Zhang @*/ 8137087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC(TS ts,PC *pc) 8147cb94ee6SHong Zhang { 8157cb94ee6SHong Zhang PetscFunctionBegin; 8165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC*),(ts,pc))); 8177cb94ee6SHong Zhang PetscFunctionReturn(0); 8187cb94ee6SHong Zhang } 8197cb94ee6SHong Zhang 820f1cd61daSJed Brown /*@ 821f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 822f1cd61daSJed Brown 823d8d19677SJose E. Roman Input Parameters: 824f1cd61daSJed Brown + ts - the time-step context 825f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 826f1cd61daSJed Brown 827fc6b9e64SJed Brown Note: 828fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 829fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 830fc6b9e64SJed Brown 831f1cd61daSJed Brown Level: beginner 832f1cd61daSJed Brown 833f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 834f1cd61daSJed Brown @*/ 8357087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 836f1cd61daSJed Brown { 837f1cd61daSJed Brown PetscFunctionBegin; 8385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt))); 839f1cd61daSJed Brown PetscFunctionReturn(0); 840f1cd61daSJed Brown } 841f1cd61daSJed Brown 842f1cd61daSJed Brown /*@ 843f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 844f1cd61daSJed Brown 845d8d19677SJose E. Roman Input Parameters: 846f1cd61daSJed Brown + ts - the time-step context 847f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 848f1cd61daSJed Brown 849f1cd61daSJed Brown Level: beginner 850f1cd61daSJed Brown 851f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 852f1cd61daSJed Brown @*/ 8537087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 854f1cd61daSJed Brown { 855f1cd61daSJed Brown PetscFunctionBegin; 8565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt))); 857f1cd61daSJed Brown PetscFunctionReturn(0); 858f1cd61daSJed Brown } 859f1cd61daSJed Brown 8602bfc04deSHong Zhang /*@ 8612bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 8622bfc04deSHong Zhang 863d8d19677SJose E. Roman Input Parameters: 8642bfc04deSHong Zhang + ts - the time-step context 8652bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 8662bfc04deSHong Zhang 8672bfc04deSHong Zhang Level: beginner 8682bfc04deSHong Zhang 869f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 8702bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 871f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8722bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 8732bfc04deSHong Zhang @*/ 8747087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft) 8752bfc04deSHong Zhang { 8762bfc04deSHong Zhang PetscFunctionBegin; 8775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft))); 8782bfc04deSHong Zhang PetscFunctionReturn(0); 8792bfc04deSHong Zhang } 880918687b8SPatrick Sanan 881918687b8SPatrick Sanan /*@ 882918687b8SPatrick Sanan TSSundialsSetUseDense - Set a flag to use a dense linear solver in SUNDIALS (serial only) 883918687b8SPatrick Sanan 884918687b8SPatrick Sanan Logically Collective 885918687b8SPatrick Sanan 886918687b8SPatrick Sanan Input Parameters: 887918687b8SPatrick Sanan + ts - the time-step context 888918687b8SPatrick Sanan - use_dense - PETSC_TRUE to use the dense solver 889918687b8SPatrick Sanan 890918687b8SPatrick Sanan Level: advanced 891918687b8SPatrick Sanan 892918687b8SPatrick Sanan .seealso: TSSUNDIALS 893918687b8SPatrick Sanan 894918687b8SPatrick Sanan @*/ 895918687b8SPatrick Sanan PetscErrorCode TSSundialsSetUseDense(TS ts,PetscBool use_dense) 896918687b8SPatrick Sanan { 897918687b8SPatrick Sanan PetscFunctionBegin; 898918687b8SPatrick Sanan PetscValidLogicalCollectiveInt(ts,use_dense,2); 8995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(ts,"TSSundialsSetUseDense_C",(TS,PetscBool),(ts,use_dense))); 900918687b8SPatrick Sanan PetscFunctionReturn(0); 901918687b8SPatrick Sanan } 902918687b8SPatrick Sanan 9037cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 9047cb94ee6SHong Zhang /*MC 90596f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 9067cb94ee6SHong Zhang 9077cb94ee6SHong Zhang Options Database: 908897c9f78SHong Zhang + -ts_sundials_type <bdf,adams> - 9097cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 9107cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 9117cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 912897c9f78SHong Zhang . -ts_sundials_linear_tolerance <tol> - 913f61b2b6aSHong Zhang . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 914918687b8SPatrick Sanan . -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps 915918687b8SPatrick Sanan - -ts_sundials_use_dense - Use a dense linear solver within CVODE (serial only) 91616016685SHong Zhang 91795452b02SPatrick Sanan Notes: 91895452b02SPatrick Sanan This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply, 919897c9f78SHong Zhang only PETSc PC options. 9207cb94ee6SHong Zhang 9217cb94ee6SHong Zhang Level: beginner 9227cb94ee6SHong Zhang 923f61b2b6aSHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(), 924a43b19c4SJed Brown TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime() 9257cb94ee6SHong Zhang 9267cb94ee6SHong Zhang M*/ 9278cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts) 9287cb94ee6SHong Zhang { 9297cb94ee6SHong Zhang TS_Sundials *cvode; 93042528757SHong Zhang PC pc; 9317cb94ee6SHong Zhang 9327cb94ee6SHong Zhang PetscFunctionBegin; 933277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Sundials; 93428adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials; 93528adb3f7SHong Zhang ts->ops->view = TSView_Sundials; 936214bc6a2SJed Brown ts->ops->setup = TSSetUp_Sundials; 937b4eba00bSSean Farley ts->ops->step = TSStep_Sundials; 938b4eba00bSSean Farley ts->ops->interpolate = TSInterpolate_Sundials; 939214bc6a2SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Sundials; 9402ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 9417cb94ee6SHong Zhang 9425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(ts,&cvode)); 943bbd56ea5SKarl Rupp 944efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 945efd4aadfSBarry Smith 9467cb94ee6SHong Zhang ts->data = (void*)cvode; 9476fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 9486fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 949f61b2b6aSHong Zhang cvode->maxl = 5; 950c4e80e11SFlorian cvode->maxord = PETSC_DEFAULT; 9517cb94ee6SHong Zhang cvode->linear_tol = .05; 952b4eba00bSSean Farley cvode->monitorstep = PETSC_TRUE; 953918687b8SPatrick Sanan cvode->use_dense = PETSC_FALSE; 9547cb94ee6SHong Zhang 9555f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)ts),&(cvode->comm_sundials))); 956f1cd61daSJed Brown 957f1cd61daSJed Brown cvode->mindt = -1.; 958f1cd61daSJed Brown cvode->maxdt = -1.; 959f1cd61daSJed Brown 9607cb94ee6SHong Zhang /* set tolerance for Sundials */ 9617cb94ee6SHong Zhang cvode->reltol = 1e-6; 9622c823083SHong Zhang cvode->abstol = 1e-6; 9637cb94ee6SHong Zhang 96442528757SHong Zhang /* set PCNONE as default pctype */ 9655f80ce2aSJacob Faibussowitsch CHKERRQ(TSSundialsGetPC_Sundials(ts,&pc)); 9665f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pc,PCNONE)); 96742528757SHong Zhang 9685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",TSSundialsSetType_Sundials)); 9695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",TSSundialsSetMaxl_Sundials)); 9705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",TSSundialsSetLinearTolerance_Sundials)); 9715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",TSSundialsSetGramSchmidtType_Sundials)); 9725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",TSSundialsSetTolerance_Sundials)); 9735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",TSSundialsSetMinTimeStep_Sundials)); 9745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",TSSundialsSetMaxTimeStep_Sundials)); 9755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",TSSundialsGetPC_Sundials)); 9765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",TSSundialsGetIterations_Sundials)); 9775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",TSSundialsMonitorInternalSteps_Sundials)); 9785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetUseDense_C",TSSundialsSetUseDense_Sundials)); 9797cb94ee6SHong Zhang PetscFunctionReturn(0); 9807cb94ee6SHong Zhang } 981