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; 207cb94ee6SHong Zhang PetscErrorCode ierr; 210679a0aeSJed Brown Mat J,P; 220679a0aeSJed Brown Vec yy = cvode->w1,yydot = cvode->ydot; 230679a0aeSJed Brown PetscReal gm = (PetscReal)_gamma; 247cb94ee6SHong Zhang PetscScalar *y_data; 257cb94ee6SHong Zhang 267cb94ee6SHong Zhang PetscFunctionBegin; 270298fd71SBarry Smith ierr = TSGetIJacobian(ts,&J,&P,NULL,NULL);CHKERRQ(ierr); 287cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(y); 297cb94ee6SHong Zhang ierr = VecPlaceArray(yy,y_data);CHKERRQ(ierr); 300679a0aeSJed Brown ierr = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */ 310679a0aeSJed Brown /* compute the shifted Jacobian (1/gm)*I + Jrest */ 3209e82e2bSBarry Smith ierr = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,J,P,PETSC_FALSE);CHKERRQ(ierr); 337cb94ee6SHong Zhang ierr = VecResetArray(yy);CHKERRQ(ierr); 340679a0aeSJed Brown ierr = MatScale(P,gm);CHKERRQ(ierr); /* turn into I-gm*Jrest, J is not used by Sundials */ 357cb94ee6SHong Zhang *jcurPtr = TRUE; 36dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 3709e82e2bSBarry Smith ierr = PCSetOperators(pc,J,P);CHKERRQ(ierr); 387cb94ee6SHong Zhang PetscFunctionReturn(0); 397cb94ee6SHong Zhang } 407cb94ee6SHong Zhang 417cb94ee6SHong Zhang /* 427cb94ee6SHong Zhang TSPSolve_Sundials - routine that we provide to Sundials that applies the preconditioner. 437cb94ee6SHong Zhang */ 444ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z, 454ac7836bSHong Zhang realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp) 467cb94ee6SHong Zhang { 477cb94ee6SHong Zhang TS ts = (TS) P_data; 487cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 49dcbc6d53SSean Farley PC pc; 507cb94ee6SHong Zhang Vec rr = cvode->w1,zz = cvode->w2; 517cb94ee6SHong Zhang PetscErrorCode ierr; 527cb94ee6SHong Zhang PetscScalar *r_data,*z_data; 537cb94ee6SHong Zhang 547cb94ee6SHong Zhang PetscFunctionBegin; 557cb94ee6SHong Zhang /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/ 567cb94ee6SHong Zhang r_data = (PetscScalar*) N_VGetArrayPointer(r); 577cb94ee6SHong Zhang z_data = (PetscScalar*) N_VGetArrayPointer(z); 587cb94ee6SHong Zhang ierr = VecPlaceArray(rr,r_data);CHKERRQ(ierr); 597cb94ee6SHong Zhang ierr = VecPlaceArray(zz,z_data);CHKERRQ(ierr); 604ac7836bSHong Zhang 617cb94ee6SHong Zhang /* Solve the Px=r and put the result in zz */ 62dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 637cb94ee6SHong Zhang ierr = PCApply(pc,rr,zz);CHKERRQ(ierr); 647cb94ee6SHong Zhang ierr = VecResetArray(rr);CHKERRQ(ierr); 657cb94ee6SHong Zhang ierr = VecResetArray(zz);CHKERRQ(ierr); 667cb94ee6SHong Zhang PetscFunctionReturn(0); 677cb94ee6SHong Zhang } 687cb94ee6SHong Zhang 697cb94ee6SHong Zhang /* 707cb94ee6SHong Zhang TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side. 717cb94ee6SHong Zhang */ 724ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx) 737cb94ee6SHong Zhang { 747cb94ee6SHong Zhang TS ts = (TS) ctx; 755fcef5e4SJed Brown DM dm; 76942e3340SBarry Smith DMTS tsdm; 775fcef5e4SJed Brown TSIFunction ifunction; 78ce94432eSBarry Smith MPI_Comm comm; 797cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 800679a0aeSJed Brown Vec yy = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot; 817cb94ee6SHong Zhang PetscScalar *y_data,*ydot_data; 827cb94ee6SHong Zhang PetscErrorCode ierr; 837cb94ee6SHong Zhang 847cb94ee6SHong Zhang PetscFunctionBegin; 85ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 865ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/ 877cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(y); 887cb94ee6SHong Zhang ydot_data = (PetscScalar*) N_VGetArrayPointer(ydot); 894ac7836bSHong Zhang ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr); 904ac7836bSHong Zhang ierr = VecPlaceArray(yyd,ydot_data);CHKERRABORT(comm,ierr); 914ac7836bSHong Zhang 925fcef5e4SJed Brown /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */ 935fcef5e4SJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 94942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 950298fd71SBarry Smith ierr = DMTSGetIFunction(dm,&ifunction,NULL);CHKERRQ(ierr); 965fcef5e4SJed Brown if (!ifunction) { 97bc0cc02bSJed Brown ierr = TSComputeRHSFunction(ts,t,yy,yyd);CHKERRQ(ierr); 98bc0cc02bSJed Brown } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */ 99bc0cc02bSJed Brown ierr = VecZeroEntries(yydot);CHKERRQ(ierr); 1000679a0aeSJed Brown ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE);CHKERRABORT(comm,ierr); 1010679a0aeSJed Brown ierr = VecScale(yyd,-1.);CHKERRQ(ierr); 102bc0cc02bSJed Brown } 1037cb94ee6SHong Zhang ierr = VecResetArray(yy);CHKERRABORT(comm,ierr); 1047cb94ee6SHong Zhang ierr = VecResetArray(yyd);CHKERRABORT(comm,ierr); 1054ac7836bSHong Zhang PetscFunctionReturn(0); 1067cb94ee6SHong Zhang } 1077cb94ee6SHong Zhang 1087cb94ee6SHong Zhang /* 109b4eba00bSSean Farley TSStep_Sundials - Calls Sundials to integrate the ODE. 1107cb94ee6SHong Zhang */ 111b4eba00bSSean Farley PetscErrorCode TSStep_Sundials(TS ts) 1127cb94ee6SHong Zhang { 1137cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 1147cb94ee6SHong Zhang PetscErrorCode ierr; 115b4eba00bSSean Farley PetscInt flag; 116be5899b3SLisandro Dalcin long int nits,lits,nsteps; 1177cb94ee6SHong Zhang realtype t,tout; 1187cb94ee6SHong Zhang PetscScalar *y_data; 1197cb94ee6SHong Zhang void *mem; 1207cb94ee6SHong Zhang 1217cb94ee6SHong Zhang PetscFunctionBegin; 12216016685SHong Zhang mem = cvode->mem; 1237cb94ee6SHong Zhang tout = ts->max_time; 1247cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr); 1257cb94ee6SHong Zhang N_VSetArrayPointer((realtype*)y_data,cvode->y); 1260298fd71SBarry Smith ierr = VecRestoreArray(ts->vec_sol,NULL);CHKERRQ(ierr); 127186e87acSLisandro Dalcin 1289687d888SLisandro Dalcin /* We would like to TSPreStage() and TSPostStage() 1299687d888SLisandro Dalcin * before each stage solve but CVode does not appear to support this. */ 130be5899b3SLisandro Dalcin if (cvode->monitorstep) 131be5899b3SLisandro Dalcin flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP); 132be5899b3SLisandro Dalcin else 133be5899b3SLisandro Dalcin flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL); 1349f94935aSHong Zhang 1359f94935aSHong Zhang if (flag) { /* display error message */ 1369f94935aSHong Zhang switch (flag) { 1379f94935aSHong Zhang case CV_ILL_INPUT: 1389f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT"); 1399f94935aSHong Zhang break; 1409f94935aSHong Zhang case CV_TOO_CLOSE: 1419f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE"); 1429f94935aSHong Zhang break; 1433c7fefeeSJed Brown case CV_TOO_MUCH_WORK: { 1449f94935aSHong Zhang PetscReal tcur; 1459f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 1469f94935aSHong Zhang ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr); 14798921bdaSJacob 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); 1483c7fefeeSJed Brown } break; 1499f94935aSHong Zhang case CV_TOO_MUCH_ACC: 1509f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC"); 1519f94935aSHong Zhang break; 1529f94935aSHong Zhang case CV_ERR_FAILURE: 1539f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE"); 1549f94935aSHong Zhang break; 1559f94935aSHong Zhang case CV_CONV_FAILURE: 1569f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE"); 1579f94935aSHong Zhang break; 1589f94935aSHong Zhang case CV_LINIT_FAIL: 1599f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL"); 1609f94935aSHong Zhang break; 1619f94935aSHong Zhang case CV_LSETUP_FAIL: 1629f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL"); 1639f94935aSHong Zhang break; 1649f94935aSHong Zhang case CV_LSOLVE_FAIL: 1659f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL"); 1669f94935aSHong Zhang break; 1679f94935aSHong Zhang case CV_RHSFUNC_FAIL: 1689f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL"); 1699f94935aSHong Zhang break; 1709f94935aSHong Zhang case CV_FIRST_RHSFUNC_ERR: 1719f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR"); 1729f94935aSHong Zhang break; 1739f94935aSHong Zhang case CV_REPTD_RHSFUNC_ERR: 1749f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR"); 1759f94935aSHong Zhang break; 1769f94935aSHong Zhang case CV_UNREC_RHSFUNC_ERR: 1779f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR"); 1789f94935aSHong Zhang break; 1799f94935aSHong Zhang case CV_RTFUNC_FAIL: 1809f94935aSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL"); 1819f94935aSHong Zhang break; 1829f94935aSHong Zhang default: 18398921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag); 1849f94935aSHong Zhang } 1859f94935aSHong Zhang } 1869f94935aSHong Zhang 187be5899b3SLisandro Dalcin /* log inner nonlinear and linear iterations */ 188be5899b3SLisandro Dalcin ierr = CVodeGetNumNonlinSolvIters(mem,&nits);CHKERRQ(ierr); 189be5899b3SLisandro Dalcin ierr = CVSpilsGetNumLinIters(mem,&lits);CHKERRQ(ierr); 190be5899b3SLisandro Dalcin ts->snes_its += nits; ts->ksp_its = lits; 191be5899b3SLisandro Dalcin 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); 196bc0cc02bSJed Brown ierr = VecCopy(cvode->update,ts->vec_sol);CHKERRQ(ierr); 197186e87acSLisandro Dalcin 198186e87acSLisandro Dalcin ts->time_step = t - ts->ptime; 199186e87acSLisandro Dalcin ts->ptime = t; 200186e87acSLisandro Dalcin 2019f94935aSHong Zhang ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr); 2022808aa04SLisandro Dalcin if (!cvode->monitorstep) ts->steps += nsteps - 1; /* TSStep() increments the step counter by one */ 203b4eba00bSSean Farley PetscFunctionReturn(0); 204b4eba00bSSean Farley } 205b4eba00bSSean Farley 206b4eba00bSSean Farley static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X) 207b4eba00bSSean Farley { 208b4eba00bSSean Farley TS_Sundials *cvode = (TS_Sundials*)ts->data; 209b4eba00bSSean Farley N_Vector y; 210b4eba00bSSean Farley PetscErrorCode ierr; 211b4eba00bSSean Farley PetscScalar *x_data; 212b4eba00bSSean Farley PetscInt glosize,locsize; 213b4eba00bSSean Farley 214b4eba00bSSean Farley PetscFunctionBegin; 215b4eba00bSSean Farley /* get the vector size */ 216b4eba00bSSean Farley ierr = VecGetSize(X,&glosize);CHKERRQ(ierr); 217b4eba00bSSean Farley ierr = VecGetLocalSize(X,&locsize);CHKERRQ(ierr); 2184d78c9e0SClaas Abert ierr = VecGetArray(X,&x_data);CHKERRQ(ierr); 219b4eba00bSSean Farley 2204d78c9e0SClaas Abert /* Initialize N_Vec y with x_data */ 221918687b8SPatrick Sanan if (cvode->use_dense) { 222918687b8SPatrick Sanan PetscMPIInt size; 223918687b8SPatrick Sanan 224918687b8SPatrick Sanan ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 225*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"TSSUNDIALS only supports a dense solve in the serial case"); 226918687b8SPatrick Sanan y = N_VMake_Serial(locsize,(realtype*)x_data); 227918687b8SPatrick Sanan } else { 2284d78c9e0SClaas Abert y = N_VMake_Parallel(cvode->comm_sundials,locsize,glosize,(realtype*)x_data); 229918687b8SPatrick Sanan } 230918687b8SPatrick Sanan 231*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!y,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Interpolated y is not allocated"); 232b4eba00bSSean Farley 233b4eba00bSSean Farley ierr = CVodeGetDky(cvode->mem,t,0,y);CHKERRQ(ierr); 234b4eba00bSSean Farley ierr = VecRestoreArray(X,&x_data);CHKERRQ(ierr); 2357cb94ee6SHong Zhang PetscFunctionReturn(0); 2367cb94ee6SHong Zhang } 2377cb94ee6SHong Zhang 238277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts) 239277b19d0SLisandro Dalcin { 240277b19d0SLisandro Dalcin TS_Sundials *cvode = (TS_Sundials*)ts->data; 241277b19d0SLisandro Dalcin PetscErrorCode ierr; 242277b19d0SLisandro Dalcin 243277b19d0SLisandro Dalcin PetscFunctionBegin; 2445c6b4a3dSJed Brown ierr = VecDestroy(&cvode->update);CHKERRQ(ierr); 2450679a0aeSJed Brown ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr); 2465c6b4a3dSJed Brown ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr); 2475c6b4a3dSJed Brown ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr); 248bbd56ea5SKarl Rupp if (cvode->mem) CVodeFree(&cvode->mem); 249277b19d0SLisandro Dalcin PetscFunctionReturn(0); 250277b19d0SLisandro Dalcin } 251277b19d0SLisandro Dalcin 2527cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts) 2537cb94ee6SHong Zhang { 2547cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2557cb94ee6SHong Zhang PetscErrorCode ierr; 2567cb94ee6SHong Zhang 2577cb94ee6SHong Zhang PetscFunctionBegin; 258277b19d0SLisandro Dalcin ierr = TSReset_Sundials(ts);CHKERRQ(ierr); 259ffc4695bSBarry Smith ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRMPI(ierr); 260277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 261bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",NULL);CHKERRQ(ierr); 262bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",NULL);CHKERRQ(ierr); 263bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",NULL);CHKERRQ(ierr); 264bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",NULL);CHKERRQ(ierr); 265bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",NULL);CHKERRQ(ierr); 266bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",NULL);CHKERRQ(ierr); 267bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",NULL);CHKERRQ(ierr); 268bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",NULL);CHKERRQ(ierr); 269bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",NULL);CHKERRQ(ierr); 270bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",NULL);CHKERRQ(ierr); 2717cb94ee6SHong Zhang PetscFunctionReturn(0); 2727cb94ee6SHong Zhang } 2737cb94ee6SHong Zhang 274214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts) 2757cb94ee6SHong Zhang { 2767cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 2777cb94ee6SHong Zhang PetscErrorCode ierr; 27816016685SHong Zhang PetscInt glosize,locsize,i,flag; 2797cb94ee6SHong Zhang PetscScalar *y_data,*parray; 28016016685SHong Zhang void *mem; 281dcbc6d53SSean Farley PC pc; 28219fd82e9SBarry Smith PCType pctype; 283ace3abfcSBarry Smith PetscBool pcnone; 2847cb94ee6SHong Zhang 2857cb94ee6SHong Zhang PetscFunctionBegin; 286*2c71b3e2SJacob 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"); 287236c45dcSShri Abhyankar 2887cb94ee6SHong Zhang /* get the vector size */ 2897cb94ee6SHong Zhang ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr); 2907cb94ee6SHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr); 2917cb94ee6SHong Zhang 2927cb94ee6SHong Zhang /* allocate the memory for N_Vec y */ 293918687b8SPatrick Sanan if (cvode->use_dense) { 294918687b8SPatrick Sanan PetscMPIInt size; 295918687b8SPatrick Sanan 296918687b8SPatrick Sanan ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 297*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"TSSUNDIALS only supports a dense solve in the serial case"); 298918687b8SPatrick Sanan cvode->y = N_VNew_Serial(locsize); 299918687b8SPatrick Sanan } else { 300a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize); 301918687b8SPatrick Sanan } 302*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!cvode->y,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cvode->y is not allocated"); 3037cb94ee6SHong Zhang 30428adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */ 3057cb94ee6SHong Zhang ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr); 3067cb94ee6SHong Zhang y_data = (PetscScalar*) N_VGetArrayPointer(cvode->y); 3077cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i]; 3080298fd71SBarry Smith ierr = VecRestoreArray(ts->vec_sol,NULL);CHKERRQ(ierr); 30942528757SHong Zhang 3107cb94ee6SHong Zhang ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr); 3110679a0aeSJed Brown ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr); 3123bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->update);CHKERRQ(ierr); 3133bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->ydot);CHKERRQ(ierr); 3147cb94ee6SHong Zhang 3157cb94ee6SHong Zhang /* 3167cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are 3177cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided 3187cb94ee6SHong Zhang by Sundials and set using VecPlaceArray(). 3197cb94ee6SHong Zhang */ 320c793f718SLisandro Dalcin ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,NULL,&cvode->w1);CHKERRQ(ierr); 321c793f718SLisandro Dalcin ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,NULL,&cvode->w2);CHKERRQ(ierr); 3223bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w1);CHKERRQ(ierr); 3233bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w2);CHKERRQ(ierr); 32416016685SHong Zhang 32516016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */ 32616016685SHong Zhang mem = CVodeCreate(cvode->cvode_type, CV_NEWTON); 327*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!mem,PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails"); 32816016685SHong Zhang cvode->mem = mem; 32916016685SHong Zhang 33016016685SHong Zhang /* Set the pointer to user-defined data */ 33116016685SHong Zhang flag = CVodeSetUserData(mem, ts); 332*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails"); 33316016685SHong Zhang 334fc6b9e64SJed Brown /* Sundials may choose to use a smaller initial step, but will never use a larger step. */ 335cdbf8f93SLisandro Dalcin flag = CVodeSetInitStep(mem,(realtype)ts->time_step); 336*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetInitStep() failed"); 337f1cd61daSJed Brown if (cvode->mindt > 0) { 338f1cd61daSJed Brown flag = CVodeSetMinStep(mem,(realtype)cvode->mindt); 3399f94935aSHong Zhang if (flag) { 340*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag == CV_MEM_NULL,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL"); 341*2c71b3e2SJacob 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"); 342ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed"); 3439f94935aSHong Zhang } 344f1cd61daSJed Brown } 345f1cd61daSJed Brown if (cvode->maxdt > 0) { 346f1cd61daSJed Brown flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt); 347*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMaxStep() failed"); 348f1cd61daSJed Brown } 349f1cd61daSJed Brown 35016016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the 351a5b23f4aSJose E. Roman * user's right hand side function in u'=f(t,u), the initial time T0, and 35216016685SHong Zhang * the initial dependent variable vector cvode->y */ 35316016685SHong Zhang flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y); 354*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag); 35516016685SHong Zhang 3569f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */ 35716016685SHong Zhang flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol); 358*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag); 35916016685SHong Zhang 360c4e80e11SFlorian /* Specify max order of BDF / ADAMS method */ 361c4e80e11SFlorian if (cvode->maxord != PETSC_DEFAULT) { 362c4e80e11SFlorian flag = CVodeSetMaxOrd(mem,cvode->maxord); 363*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxOrd() fails, flag %d",flag); 364c4e80e11SFlorian } 365c4e80e11SFlorian 3669f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */ 3679f94935aSHong Zhang flag = CVodeSetMaxNumSteps(mem,ts->max_steps); 368*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxNumSteps() fails, flag %d",flag); 36916016685SHong Zhang 370918687b8SPatrick Sanan if (cvode->use_dense) { 371918687b8SPatrick Sanan /* call CVDense to use a dense linear solver. */ 372918687b8SPatrick Sanan PetscMPIInt size; 373918687b8SPatrick Sanan 374918687b8SPatrick Sanan ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 375*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"TSSUNDIALS only supports a dense solve in the serial case"); 376918687b8SPatrick Sanan flag = CVDense(mem,locsize); 377*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVDense() fails, flag %d",flag); 378918687b8SPatrick Sanan } else { 37916016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */ 38016016685SHong Zhang /* setup the ode integrator with the given preconditioner */ 381dcbc6d53SSean Farley ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 382dcbc6d53SSean Farley ierr = PCGetType(pc,&pctype);CHKERRQ(ierr); 383251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr); 38416016685SHong Zhang if (pcnone) { 38516016685SHong Zhang flag = CVSpgmr(mem,PREC_NONE,0); 386*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 38716016685SHong Zhang } else { 388f61b2b6aSHong Zhang flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl); 389*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag); 39016016685SHong Zhang 39116016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve, 39216016685SHong Zhang and the pointer to the user-defined block data */ 39316016685SHong Zhang flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials); 394*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag); 39516016685SHong Zhang } 396918687b8SPatrick Sanan } 3977cb94ee6SHong Zhang PetscFunctionReturn(0); 3987cb94ee6SHong Zhang } 3997cb94ee6SHong Zhang 4006fadb2cdSHong Zhang /* type of CVODE linear multistep method */ 401c793f718SLisandro Dalcin const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",NULL}; 4026fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */ 403c793f718SLisandro Dalcin const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",NULL}; 404a04cf4d8SBarry Smith 4054416b707SBarry Smith PetscErrorCode TSSetFromOptions_Sundials(PetscOptionItems *PetscOptionsObject,TS ts) 4067cb94ee6SHong Zhang { 4077cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4087cb94ee6SHong Zhang PetscErrorCode ierr; 4097cb94ee6SHong Zhang int indx; 410ace3abfcSBarry Smith PetscBool flag; 4117bda3a07SJed Brown PC pc; 4127cb94ee6SHong Zhang 4137cb94ee6SHong Zhang PetscFunctionBegin; 414e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SUNDIALS ODE solver options");CHKERRQ(ierr); 4156fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr); 4167cb94ee6SHong Zhang if (flag) { 4176fadb2cdSHong Zhang ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr); 4187cb94ee6SHong Zhang } 4196fadb2cdSHong Zhang ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr); 4207cb94ee6SHong Zhang if (flag) { 4217cb94ee6SHong Zhang ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr); 4227cb94ee6SHong Zhang } 4230298fd71SBarry Smith ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL);CHKERRQ(ierr); 4240298fd71SBarry Smith ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL);CHKERRQ(ierr); 4250298fd71SBarry Smith ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL);CHKERRQ(ierr); 4260298fd71SBarry Smith ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL);CHKERRQ(ierr); 42794ae4db5SBarry Smith ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,NULL);CHKERRQ(ierr); 428c4e80e11SFlorian ierr = PetscOptionsInt("-ts_sundials_maxord","Max Order for BDF/Adams method","TSSundialsSetMaxOrd",cvode->maxord,&cvode->maxord,NULL);CHKERRQ(ierr); 42994ae4db5SBarry Smith ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,NULL);CHKERRQ(ierr); 430be5899b3SLisandro Dalcin ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internal steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL);CHKERRQ(ierr); 431918687b8SPatrick Sanan ierr = PetscOptionsBool("-ts_sundials_use_dense","Use dense internal solver in SUNDIALS (serial only)","TSSundialsSetUseDense",cvode->use_dense,&cvode->use_dense,NULL);CHKERRQ(ierr); 4327cb94ee6SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4337bda3a07SJed Brown ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 4347bda3a07SJed Brown ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 4357cb94ee6SHong Zhang PetscFunctionReturn(0); 4367cb94ee6SHong Zhang } 4377cb94ee6SHong Zhang 4387cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer) 4397cb94ee6SHong Zhang { 4407cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 4417cb94ee6SHong Zhang PetscErrorCode ierr; 4427cb94ee6SHong Zhang char *type; 4437cb94ee6SHong Zhang char atype[] = "Adams"; 4447cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula"; 445ace3abfcSBarry Smith PetscBool iascii,isstring; 4462c823083SHong Zhang long int nsteps,its,nfevals,nlinsetups,nfails,itmp; 4472c823083SHong Zhang PetscInt qlast,qcur; 4485d47aee6SHong Zhang PetscReal hinused,hlast,hcur,tcur,tolsfac; 44942528757SHong Zhang PC pc; 4507cb94ee6SHong Zhang 4517cb94ee6SHong Zhang PetscFunctionBegin; 452bbd56ea5SKarl Rupp if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype; 453bbd56ea5SKarl Rupp else type = btype; 4547cb94ee6SHong Zhang 455251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 456251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 4577cb94ee6SHong Zhang if (iascii) { 458918687b8SPatrick Sanan ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrator does not use SNES!\n");CHKERRQ(ierr); 459918687b8SPatrick Sanan ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrator type %s\n",type);CHKERRQ(ierr); 460c4e80e11SFlorian ierr = PetscViewerASCIIPrintf(viewer,"Sundials maxord %D\n",cvode->maxord);CHKERRQ(ierr); 461c4e80e11SFlorian ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",(double)cvode->abstol,(double)cvode->reltol);CHKERRQ(ierr); 462918687b8SPatrick Sanan if (cvode->use_dense) { 463918687b8SPatrick Sanan ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrator using a dense linear solve\n");CHKERRQ(ierr); 464918687b8SPatrick Sanan } else { 465c4e80e11SFlorian ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",(double)cvode->linear_tol);CHKERRQ(ierr); 466f61b2b6aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr); 4677cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) { 4687cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4697cb94ee6SHong Zhang } else { 4707cb94ee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr); 4717cb94ee6SHong Zhang } 472918687b8SPatrick Sanan } 473c4e80e11SFlorian if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",(double)cvode->mindt);CHKERRQ(ierr);} 474c4e80e11SFlorian if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",(double)cvode->maxdt);CHKERRQ(ierr);} 4752c823083SHong Zhang 4765d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */ 4775d47aee6SHong Zhang ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr); 4785d47aee6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr); 4792c823083SHong Zhang ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals, 4802c823083SHong Zhang &nlinsetups,&nfails,&qlast,&qcur, 4812c823083SHong Zhang &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr); 4822c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr); 4832c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr); 4842c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr); 4852c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr); 4862c823083SHong Zhang 4872c823083SHong Zhang ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr); 4882c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr); 4892c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr); 490918687b8SPatrick Sanan if (!cvode->use_dense) { 4912c823083SHong Zhang ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */ 4922c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr); 4932c823083SHong Zhang ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr); 4942c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr); 49542528757SHong Zhang 49642528757SHong Zhang ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr); 49742528757SHong Zhang ierr = PCView(pc,viewer);CHKERRQ(ierr); 4982c823083SHong Zhang ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr); 4992c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr); 5002c823083SHong Zhang ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr); 5012c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr); 502918687b8SPatrick Sanan } 5032c823083SHong Zhang ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr); 5042c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr); 5052c823083SHong Zhang ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr); 5062c823083SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr); 5077cb94ee6SHong Zhang } else if (isstring) { 5087cb94ee6SHong Zhang ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr); 5097cb94ee6SHong Zhang } 5107cb94ee6SHong Zhang PetscFunctionReturn(0); 5117cb94ee6SHong Zhang } 5127cb94ee6SHong Zhang 5137cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/ 5147087cfbeSBarry Smith PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type) 5157cb94ee6SHong Zhang { 5167cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5177cb94ee6SHong Zhang 5187cb94ee6SHong Zhang PetscFunctionBegin; 5197cb94ee6SHong Zhang cvode->cvode_type = type; 5207cb94ee6SHong Zhang PetscFunctionReturn(0); 5217cb94ee6SHong Zhang } 5227cb94ee6SHong Zhang 523f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl) 5247cb94ee6SHong Zhang { 5257cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5267cb94ee6SHong Zhang 5277cb94ee6SHong Zhang PetscFunctionBegin; 528f61b2b6aSHong Zhang cvode->maxl = maxl; 5297cb94ee6SHong Zhang PetscFunctionReturn(0); 5307cb94ee6SHong Zhang } 5317cb94ee6SHong Zhang 532590ae059SPatrick Sanan PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,PetscReal tol) 5337cb94ee6SHong Zhang { 5347cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5357cb94ee6SHong Zhang 5367cb94ee6SHong Zhang PetscFunctionBegin; 5377cb94ee6SHong Zhang cvode->linear_tol = tol; 5387cb94ee6SHong Zhang PetscFunctionReturn(0); 5397cb94ee6SHong Zhang } 5407cb94ee6SHong Zhang 5417087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type) 5427cb94ee6SHong Zhang { 5437cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5447cb94ee6SHong Zhang 5457cb94ee6SHong Zhang PetscFunctionBegin; 5467cb94ee6SHong Zhang cvode->gtype = type; 5477cb94ee6SHong Zhang PetscFunctionReturn(0); 5487cb94ee6SHong Zhang } 5497cb94ee6SHong Zhang 550590ae059SPatrick Sanan PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,PetscReal aabs,PetscReal rel) 5517cb94ee6SHong Zhang { 5527cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 5537cb94ee6SHong Zhang 5547cb94ee6SHong Zhang PetscFunctionBegin; 5557cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs; 5567cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel; 5577cb94ee6SHong Zhang PetscFunctionReturn(0); 5587cb94ee6SHong Zhang } 5597cb94ee6SHong Zhang 5607087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt) 561f1cd61daSJed Brown { 562f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 563f1cd61daSJed Brown 564f1cd61daSJed Brown PetscFunctionBegin; 565f1cd61daSJed Brown cvode->mindt = mindt; 566f1cd61daSJed Brown PetscFunctionReturn(0); 567f1cd61daSJed Brown } 568f1cd61daSJed Brown 5697087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt) 570f1cd61daSJed Brown { 571f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials*)ts->data; 572f1cd61daSJed Brown 573f1cd61daSJed Brown PetscFunctionBegin; 574f1cd61daSJed Brown cvode->maxdt = maxdt; 575f1cd61daSJed Brown PetscFunctionReturn(0); 576f1cd61daSJed Brown } 577918687b8SPatrick Sanan 578918687b8SPatrick Sanan PetscErrorCode TSSundialsSetUseDense_Sundials(TS ts,PetscBool use_dense) 579918687b8SPatrick Sanan { 580918687b8SPatrick Sanan TS_Sundials *cvode = (TS_Sundials*)ts->data; 581918687b8SPatrick Sanan 582918687b8SPatrick Sanan PetscFunctionBegin; 583918687b8SPatrick Sanan cvode->use_dense = use_dense; 584918687b8SPatrick Sanan PetscFunctionReturn(0); 585918687b8SPatrick Sanan } 586918687b8SPatrick Sanan 5877087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc) 5887cb94ee6SHong Zhang { 589dcbc6d53SSean Farley SNES snes; 590dcbc6d53SSean Farley KSP ksp; 591dcbc6d53SSean Farley PetscErrorCode ierr; 5927cb94ee6SHong Zhang 5937cb94ee6SHong Zhang PetscFunctionBegin; 594dcbc6d53SSean Farley ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 595dcbc6d53SSean Farley ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 596dcbc6d53SSean Farley ierr = KSPGetPC(ksp,pc);CHKERRQ(ierr); 5977cb94ee6SHong Zhang PetscFunctionReturn(0); 5987cb94ee6SHong Zhang } 5997cb94ee6SHong Zhang 6007087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin) 6017cb94ee6SHong Zhang { 6027cb94ee6SHong Zhang PetscFunctionBegin; 6035ef26d82SJed Brown if (nonlin) *nonlin = ts->snes_its; 6045ef26d82SJed Brown if (lin) *lin = ts->ksp_its; 6057cb94ee6SHong Zhang PetscFunctionReturn(0); 6067cb94ee6SHong Zhang } 6077cb94ee6SHong Zhang 6087087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s) 6092bfc04deSHong Zhang { 6102bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials*)ts->data; 6112bfc04deSHong Zhang 6122bfc04deSHong Zhang PetscFunctionBegin; 6132bfc04deSHong Zhang cvode->monitorstep = s; 6142bfc04deSHong Zhang PetscFunctionReturn(0); 6152bfc04deSHong Zhang } 6167cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 6177cb94ee6SHong Zhang 6187cb94ee6SHong Zhang /*@C 6197cb94ee6SHong Zhang TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials. 6207cb94ee6SHong Zhang 6217cb94ee6SHong Zhang Not Collective 6227cb94ee6SHong Zhang 62397bb3fdcSJose E. Roman Input Parameter: 6247cb94ee6SHong Zhang . ts - the time-step context 6257cb94ee6SHong Zhang 6267cb94ee6SHong Zhang Output Parameters: 6277cb94ee6SHong Zhang + nonlin - number of nonlinear iterations 6287cb94ee6SHong Zhang - lin - number of linear iterations 6297cb94ee6SHong Zhang 6307cb94ee6SHong Zhang Level: advanced 6317cb94ee6SHong Zhang 6327cb94ee6SHong Zhang Notes: 6337cb94ee6SHong Zhang These return the number since the creation of the TS object 6347cb94ee6SHong Zhang 635f61b2b6aSHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetMaxl(), 6367cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 637f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 638a43b19c4SJed Brown TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime() 6397cb94ee6SHong Zhang 6407cb94ee6SHong Zhang @*/ 6417087cfbeSBarry Smith PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin) 6427cb94ee6SHong Zhang { 6434ac538c5SBarry Smith PetscErrorCode ierr; 6447cb94ee6SHong Zhang 6457cb94ee6SHong Zhang PetscFunctionBegin; 6464ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr); 6477cb94ee6SHong Zhang PetscFunctionReturn(0); 6487cb94ee6SHong Zhang } 6497cb94ee6SHong Zhang 6507cb94ee6SHong Zhang /*@ 6517cb94ee6SHong Zhang TSSundialsSetType - Sets the method that Sundials will use for integration. 6527cb94ee6SHong Zhang 6533f9fe445SBarry Smith Logically Collective on TS 6547cb94ee6SHong Zhang 6558f7d6fe5SPatrick Sanan Input Parameters: 6567cb94ee6SHong Zhang + ts - the time-step context 6577cb94ee6SHong Zhang - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF 6587cb94ee6SHong Zhang 6597cb94ee6SHong Zhang Level: intermediate 6607cb94ee6SHong Zhang 661f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetMaxl(), 6627cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 663f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 6647cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 665a43b19c4SJed Brown TSSetExactFinalTime() 6667cb94ee6SHong Zhang @*/ 6677087cfbeSBarry Smith PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type) 6687cb94ee6SHong Zhang { 6694ac538c5SBarry Smith PetscErrorCode ierr; 6707cb94ee6SHong Zhang 6717cb94ee6SHong Zhang PetscFunctionBegin; 6724ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr); 6737cb94ee6SHong Zhang PetscFunctionReturn(0); 6747cb94ee6SHong Zhang } 6757cb94ee6SHong Zhang 6767cb94ee6SHong Zhang /*@ 677c4e80e11SFlorian TSSundialsSetMaxord - Sets the maximum order for BDF/Adams method used by SUNDIALS. 678c4e80e11SFlorian 679c4e80e11SFlorian Logically Collective on TS 680c4e80e11SFlorian 6818f7d6fe5SPatrick Sanan Input Parameters: 682c4e80e11SFlorian + ts - the time-step context 683c4e80e11SFlorian - maxord - maximum order of BDF / Adams method 684c4e80e11SFlorian 685c4e80e11SFlorian Level: advanced 686c4e80e11SFlorian 687c4e80e11SFlorian .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 688c4e80e11SFlorian TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 689c4e80e11SFlorian TSSundialsGetIterations(), TSSundialsSetType(), 690c4e80e11SFlorian TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 691c4e80e11SFlorian TSSetExactFinalTime() 692c4e80e11SFlorian 693c4e80e11SFlorian @*/ 694c4e80e11SFlorian PetscErrorCode TSSundialsSetMaxord(TS ts,PetscInt maxord) 695c4e80e11SFlorian { 696c4e80e11SFlorian PetscErrorCode ierr; 697c4e80e11SFlorian 698c4e80e11SFlorian PetscFunctionBegin; 699c4e80e11SFlorian PetscValidLogicalCollectiveInt(ts,maxord,2); 700c4e80e11SFlorian ierr = PetscTryMethod(ts,"TSSundialsSetMaxOrd_C",(TS,PetscInt),(ts,maxord));CHKERRQ(ierr); 701c4e80e11SFlorian PetscFunctionReturn(0); 702c4e80e11SFlorian } 703c4e80e11SFlorian 704c4e80e11SFlorian /*@ 705f61b2b6aSHong Zhang TSSundialsSetMaxl - Sets the dimension of the Krylov space used by 7067cb94ee6SHong Zhang GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so 707f61b2b6aSHong Zhang this is the maximum number of GMRES steps that will be used. 7087cb94ee6SHong Zhang 7093f9fe445SBarry Smith Logically Collective on TS 7107cb94ee6SHong Zhang 7118f7d6fe5SPatrick Sanan Input Parameters: 7127cb94ee6SHong Zhang + ts - the time-step context 713f61b2b6aSHong Zhang - maxl - number of direction vectors (the dimension of Krylov subspace). 7147cb94ee6SHong Zhang 7157cb94ee6SHong Zhang Level: advanced 7167cb94ee6SHong Zhang 7177cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 7187cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 719f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7207cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 721a43b19c4SJed Brown TSSetExactFinalTime() 7227cb94ee6SHong Zhang 7237cb94ee6SHong Zhang @*/ 724f61b2b6aSHong Zhang PetscErrorCode TSSundialsSetMaxl(TS ts,PetscInt maxl) 7257cb94ee6SHong Zhang { 7264ac538c5SBarry Smith PetscErrorCode ierr; 7277cb94ee6SHong Zhang 7287cb94ee6SHong Zhang PetscFunctionBegin; 729f61b2b6aSHong Zhang PetscValidLogicalCollectiveInt(ts,maxl,2); 730f61b2b6aSHong Zhang ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));CHKERRQ(ierr); 7317cb94ee6SHong Zhang PetscFunctionReturn(0); 7327cb94ee6SHong Zhang } 7337cb94ee6SHong Zhang 7347cb94ee6SHong Zhang /*@ 7357cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear 7367cb94ee6SHong Zhang system by SUNDIALS. 7377cb94ee6SHong Zhang 7383f9fe445SBarry Smith Logically Collective on TS 7397cb94ee6SHong Zhang 7408f7d6fe5SPatrick Sanan Input Parameters: 7417cb94ee6SHong Zhang + ts - the time-step context 7427cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is 7437cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default. 7447cb94ee6SHong Zhang 7457cb94ee6SHong Zhang Level: advanced 7467cb94ee6SHong Zhang 747f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7487cb94ee6SHong Zhang TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 749f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7507cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 751a43b19c4SJed Brown TSSetExactFinalTime() 7527cb94ee6SHong Zhang 7537cb94ee6SHong Zhang @*/ 754590ae059SPatrick Sanan PetscErrorCode TSSundialsSetLinearTolerance(TS ts,PetscReal tol) 7557cb94ee6SHong Zhang { 7564ac538c5SBarry Smith PetscErrorCode ierr; 7577cb94ee6SHong Zhang 7587cb94ee6SHong Zhang PetscFunctionBegin; 759c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,tol,2); 760590ae059SPatrick Sanan ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,PetscReal),(ts,tol));CHKERRQ(ierr); 7617cb94ee6SHong Zhang PetscFunctionReturn(0); 7627cb94ee6SHong Zhang } 7637cb94ee6SHong Zhang 7647cb94ee6SHong Zhang /*@ 7657cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used 7667cb94ee6SHong Zhang in GMRES method by SUNDIALS linear solver. 7677cb94ee6SHong Zhang 7683f9fe445SBarry Smith Logically Collective on TS 7697cb94ee6SHong Zhang 7708f7d6fe5SPatrick Sanan Input Parameters: 7717cb94ee6SHong Zhang + ts - the time-step context 7727cb94ee6SHong Zhang - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS 7737cb94ee6SHong Zhang 7747cb94ee6SHong Zhang Level: advanced 7757cb94ee6SHong Zhang 776f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 7777cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), 778f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 7797cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 780a43b19c4SJed Brown TSSetExactFinalTime() 7817cb94ee6SHong Zhang 7827cb94ee6SHong Zhang @*/ 7837087cfbeSBarry Smith PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type) 7847cb94ee6SHong Zhang { 7854ac538c5SBarry Smith PetscErrorCode ierr; 7867cb94ee6SHong Zhang 7877cb94ee6SHong Zhang PetscFunctionBegin; 7884ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr); 7897cb94ee6SHong Zhang PetscFunctionReturn(0); 7907cb94ee6SHong Zhang } 7917cb94ee6SHong Zhang 7927cb94ee6SHong Zhang /*@ 7937cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 7947cb94ee6SHong Zhang Sundials for error control. 7957cb94ee6SHong Zhang 7963f9fe445SBarry Smith Logically Collective on TS 7977cb94ee6SHong Zhang 7988f7d6fe5SPatrick Sanan Input Parameters: 7997cb94ee6SHong Zhang + ts - the time-step context 8007cb94ee6SHong Zhang . aabs - the absolute tolerance 8017cb94ee6SHong Zhang - rel - the relative tolerance 8027cb94ee6SHong Zhang 8037cb94ee6SHong Zhang See the Cvode/Sundials users manual for exact details on these parameters. Essentially 8047cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep. 8057cb94ee6SHong Zhang 8067cb94ee6SHong Zhang Level: intermediate 8077cb94ee6SHong Zhang 808f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(), 8097cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 810f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8117cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(), 812a43b19c4SJed Brown TSSetExactFinalTime() 8137cb94ee6SHong Zhang 8147cb94ee6SHong Zhang @*/ 815590ae059SPatrick Sanan PetscErrorCode TSSundialsSetTolerance(TS ts,PetscReal aabs,PetscReal rel) 8167cb94ee6SHong Zhang { 8174ac538c5SBarry Smith PetscErrorCode ierr; 8187cb94ee6SHong Zhang 8197cb94ee6SHong Zhang PetscFunctionBegin; 820590ae059SPatrick Sanan ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,PetscReal,PetscReal),(ts,aabs,rel));CHKERRQ(ierr); 8217cb94ee6SHong Zhang PetscFunctionReturn(0); 8227cb94ee6SHong Zhang } 8237cb94ee6SHong Zhang 8247cb94ee6SHong Zhang /*@ 8257cb94ee6SHong Zhang TSSundialsGetPC - Extract the PC context from a time-step context for Sundials. 8267cb94ee6SHong Zhang 8277cb94ee6SHong Zhang Input Parameter: 8287cb94ee6SHong Zhang . ts - the time-step context 8297cb94ee6SHong Zhang 8307cb94ee6SHong Zhang Output Parameter: 8317cb94ee6SHong Zhang . pc - the preconditioner context 8327cb94ee6SHong Zhang 8337cb94ee6SHong Zhang Level: advanced 8347cb94ee6SHong Zhang 835f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 8367cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 837f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 8387cb94ee6SHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance() 8397cb94ee6SHong Zhang @*/ 8407087cfbeSBarry Smith PetscErrorCode TSSundialsGetPC(TS ts,PC *pc) 8417cb94ee6SHong Zhang { 8424ac538c5SBarry Smith PetscErrorCode ierr; 8437cb94ee6SHong Zhang 8447cb94ee6SHong Zhang PetscFunctionBegin; 8454ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC*),(ts,pc));CHKERRQ(ierr); 8467cb94ee6SHong Zhang PetscFunctionReturn(0); 8477cb94ee6SHong Zhang } 8487cb94ee6SHong Zhang 849f1cd61daSJed Brown /*@ 850f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller. 851f1cd61daSJed Brown 852d8d19677SJose E. Roman Input Parameters: 853f1cd61daSJed Brown + ts - the time-step context 854f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate 855f1cd61daSJed Brown 856fc6b9e64SJed Brown Note: 857fc6b9e64SJed Brown Sundials will error if it is not possible to keep the estimated truncation error below 858fc6b9e64SJed Brown the tolerance set with TSSundialsSetTolerance() without going below this step size. 859fc6b9e64SJed Brown 860f1cd61daSJed Brown Level: beginner 861f1cd61daSJed Brown 862f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 863f1cd61daSJed Brown @*/ 8647087cfbeSBarry Smith PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt) 865f1cd61daSJed Brown { 8664ac538c5SBarry Smith PetscErrorCode ierr; 867f1cd61daSJed Brown 868f1cd61daSJed Brown PetscFunctionBegin; 8694ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr); 870f1cd61daSJed Brown PetscFunctionReturn(0); 871f1cd61daSJed Brown } 872f1cd61daSJed Brown 873f1cd61daSJed Brown /*@ 874f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller. 875f1cd61daSJed Brown 876d8d19677SJose E. Roman Input Parameters: 877f1cd61daSJed Brown + ts - the time-step context 878f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate 879f1cd61daSJed Brown 880f1cd61daSJed Brown Level: beginner 881f1cd61daSJed Brown 882f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(), 883f1cd61daSJed Brown @*/ 8847087cfbeSBarry Smith PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt) 885f1cd61daSJed Brown { 8864ac538c5SBarry Smith PetscErrorCode ierr; 887f1cd61daSJed Brown 888f1cd61daSJed Brown PetscFunctionBegin; 8894ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 890f1cd61daSJed Brown PetscFunctionReturn(0); 891f1cd61daSJed Brown } 892f1cd61daSJed Brown 8932bfc04deSHong Zhang /*@ 8942bfc04deSHong Zhang TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false). 8952bfc04deSHong Zhang 896d8d19677SJose E. Roman Input Parameters: 8972bfc04deSHong Zhang + ts - the time-step context 8982bfc04deSHong Zhang - ft - PETSC_TRUE if monitor, else PETSC_FALSE 8992bfc04deSHong Zhang 9002bfc04deSHong Zhang Level: beginner 9012bfc04deSHong Zhang 902f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(), 9032bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), 904f61b2b6aSHong Zhang TSSundialsGetIterations(), TSSundialsSetType(), 9052bfc04deSHong Zhang TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 9062bfc04deSHong Zhang @*/ 9077087cfbeSBarry Smith PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft) 9082bfc04deSHong Zhang { 9094ac538c5SBarry Smith PetscErrorCode ierr; 9102bfc04deSHong Zhang 9112bfc04deSHong Zhang PetscFunctionBegin; 9124ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr); 9132bfc04deSHong Zhang PetscFunctionReturn(0); 9142bfc04deSHong Zhang } 915918687b8SPatrick Sanan 916918687b8SPatrick Sanan /*@ 917918687b8SPatrick Sanan TSSundialsSetUseDense - Set a flag to use a dense linear solver in SUNDIALS (serial only) 918918687b8SPatrick Sanan 919918687b8SPatrick Sanan Logically Collective 920918687b8SPatrick Sanan 921918687b8SPatrick Sanan Input Parameters: 922918687b8SPatrick Sanan + ts - the time-step context 923918687b8SPatrick Sanan - use_dense - PETSC_TRUE to use the dense solver 924918687b8SPatrick Sanan 925918687b8SPatrick Sanan Level: advanced 926918687b8SPatrick Sanan 927918687b8SPatrick Sanan .seealso: TSSUNDIALS 928918687b8SPatrick Sanan 929918687b8SPatrick Sanan @*/ 930918687b8SPatrick Sanan PetscErrorCode TSSundialsSetUseDense(TS ts,PetscBool use_dense) 931918687b8SPatrick Sanan { 932918687b8SPatrick Sanan PetscErrorCode ierr; 933918687b8SPatrick Sanan 934918687b8SPatrick Sanan PetscFunctionBegin; 935918687b8SPatrick Sanan PetscValidLogicalCollectiveInt(ts,use_dense,2); 936918687b8SPatrick Sanan ierr = PetscTryMethod(ts,"TSSundialsSetUseDense_C",(TS,PetscBool),(ts,use_dense));CHKERRQ(ierr); 937918687b8SPatrick Sanan PetscFunctionReturn(0); 938918687b8SPatrick Sanan } 939918687b8SPatrick Sanan 9407cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/ 9417cb94ee6SHong Zhang /*MC 94296f5712cSJed Brown TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS) 9437cb94ee6SHong Zhang 9447cb94ee6SHong Zhang Options Database: 945897c9f78SHong Zhang + -ts_sundials_type <bdf,adams> - 9467cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES 9477cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence 9487cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence 949897c9f78SHong Zhang . -ts_sundials_linear_tolerance <tol> - 950f61b2b6aSHong Zhang . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace 951918687b8SPatrick Sanan . -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps 952918687b8SPatrick Sanan - -ts_sundials_use_dense - Use a dense linear solver within CVODE (serial only) 95316016685SHong Zhang 95495452b02SPatrick Sanan Notes: 95595452b02SPatrick Sanan This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply, 956897c9f78SHong Zhang only PETSc PC options. 9577cb94ee6SHong Zhang 9587cb94ee6SHong Zhang Level: beginner 9597cb94ee6SHong Zhang 960f61b2b6aSHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(), 961a43b19c4SJed Brown TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime() 9627cb94ee6SHong Zhang 9637cb94ee6SHong Zhang M*/ 9648cc058d9SJed Brown PETSC_EXTERN 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; 9782ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 9797cb94ee6SHong Zhang 980b00a9115SJed Brown ierr = PetscNewLog(ts,&cvode);CHKERRQ(ierr); 981bbd56ea5SKarl Rupp 982efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 983efd4aadfSBarry Smith 9847cb94ee6SHong Zhang ts->data = (void*)cvode; 9856fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF; 9866fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS; 987f61b2b6aSHong Zhang cvode->maxl = 5; 988c4e80e11SFlorian cvode->maxord = PETSC_DEFAULT; 9897cb94ee6SHong Zhang cvode->linear_tol = .05; 990b4eba00bSSean Farley cvode->monitorstep = PETSC_TRUE; 991918687b8SPatrick Sanan cvode->use_dense = PETSC_FALSE; 9927cb94ee6SHong Zhang 993ffc4695bSBarry Smith ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)ts),&(cvode->comm_sundials));CHKERRMPI(ierr); 994f1cd61daSJed Brown 995f1cd61daSJed Brown cvode->mindt = -1.; 996f1cd61daSJed Brown cvode->maxdt = -1.; 997f1cd61daSJed Brown 9987cb94ee6SHong Zhang /* set tolerance for Sundials */ 9997cb94ee6SHong Zhang cvode->reltol = 1e-6; 10002c823083SHong Zhang cvode->abstol = 1e-6; 10017cb94ee6SHong Zhang 100242528757SHong Zhang /* set PCNONE as default pctype */ 100342528757SHong Zhang ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr); 100442528757SHong Zhang ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 100542528757SHong Zhang 1006bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",TSSundialsSetType_Sundials);CHKERRQ(ierr); 1007bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",TSSundialsSetMaxl_Sundials);CHKERRQ(ierr); 1008bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr); 1009bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr); 1010bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",TSSundialsSetTolerance_Sundials);CHKERRQ(ierr); 1011bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr); 1012bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr); 1013bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",TSSundialsGetPC_Sundials);CHKERRQ(ierr); 1014bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",TSSundialsGetIterations_Sundials);CHKERRQ(ierr); 1015bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr); 1016918687b8SPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetUseDense_C",TSSundialsSetUseDense_Sundials);CHKERRQ(ierr); 10177cb94ee6SHong Zhang PetscFunctionReturn(0); 10187cb94ee6SHong Zhang } 1019