82 virtual void ImplicitSolve(
const double dt,
const Vector &u,
Vector &k);
86 void SundialsSolve(
const double dt,
Vector &b);
89 void SetParameters(
const Vector &u);
91 virtual ~ConductionOperator();
107 ConductionOperator *oper;
110 SundialsJacSolver() : oper(NULL) { }
112 int InitSystem(
void *sundials_mem);
113 int SetupSystem(
void *sundials_mem,
int conv_fail,
114 const Vector &y_pred,
const Vector &f_pred,
int &jac_cur,
116 int SolveSystem(
void *sundials_mem,
Vector &b,
const Vector &weight,
118 int FreeSystem(
void *sundials_mem);
123 int main(
int argc,
char *argv[])
127 MPI_Init(&argc, &argv);
128 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
129 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
132 const char *mesh_file =
"../../data/star.mesh";
133 int ser_ref_levels = 2;
134 int par_ref_levels = 1;
136 int ode_solver_type = 11;
137 double t_final = 0.5;
139 double alpha = 1.0e-2;
141 bool visualization =
true;
146 const double reltol = 1e-4, abstol = 1e-4;
149 cout.precision(precision);
152 args.
AddOption(&mesh_file,
"-m",
"--mesh",
153 "Mesh file to use.");
154 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
155 "Number of times to refine the mesh uniformly in serial.");
156 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
157 "Number of times to refine the mesh uniformly in parallel.");
159 "Order (degree) of the finite elements.");
160 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
162 "\t 1/11 - CVODE (explicit/implicit),\n"
163 "\t 2/12 - ARKODE (default explicit/implicit),\n"
164 "\t 3 - ARKODE (Fehlberg-6-4-5)\n"
165 "\t 4 - Forward Euler, 5 - RK2, 6 - RK3 SSP, 7 - RK4,\n"
166 "\t 8 - Backward Euler, 9 - SDIRK23, 10 - SDIRK33.");
167 args.
AddOption(&t_final,
"-tf",
"--t-final",
168 "Final time; start time is 0.");
169 args.
AddOption(&dt,
"-dt",
"--time-step",
172 "Alpha coefficient.");
174 "Kappa coefficient offset.");
175 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
176 "--no-visualization",
177 "Enable or disable GLVis visualization.");
178 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
179 "--no-visit-datafiles",
180 "Save data files for VisIt (visit.llnl.gov) visualization.");
181 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
182 "Visualize every n-th timestep.");
199 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
208 SundialsJacSolver sun_solver;
209 switch (ode_solver_type)
213 cvode =
new CVODESolver(MPI_COMM_WORLD, CV_ADAMS, CV_FUNCTIONAL);
216 ode_solver = cvode;
break;
218 cvode =
new CVODESolver(MPI_COMM_WORLD, CV_BDF, CV_NEWTON);
222 ode_solver = cvode;
break;
225 arkode =
new ARKODESolver(MPI_COMM_WORLD, ARKODESolver::EXPLICIT);
228 if (ode_solver_type == 3) { arkode->
SetERKTableNum(FEHLBERG_13_7_8); }
229 ode_solver = arkode;
break;
231 arkode =
new ARKODESolver(MPI_COMM_WORLD, ARKODESolver::IMPLICIT);
235 ode_solver = arkode;
break;
238 case 5: ode_solver =
new RK2Solver(0.5);
break;
240 case 7: ode_solver =
new RK4Solver;
break;
246 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
259 for (
int lev = 0; lev < ser_ref_levels; lev++)
269 for (
int lev = 0; lev < par_ref_levels; lev++)
282 cout <<
"Number of temperature unknowns: " << fe_size << endl;
295 ConductionOperator oper(fespace, alpha, kappa, u);
299 ostringstream mesh_name, sol_name;
300 mesh_name <<
"ex16-mesh." << setfill(
'0') << setw(6) << myid;
301 sol_name <<
"ex16-init." << setfill(
'0') << setw(6) << myid;
302 ofstream omesh(mesh_name.str().c_str());
303 omesh.precision(precision);
305 ofstream osol(sol_name.str().c_str());
306 osol.precision(precision);
322 char vishost[] =
"localhost";
324 sout.
open(vishost, visport);
325 sout <<
"parallel " << num_procs <<
" " << myid << endl;
326 int good = sout.good(), all_good;
327 MPI_Allreduce(&good, &all_good, 1, MPI_INT, MPI_MIN, pmesh->
GetComm());
331 visualization =
false;
334 cout <<
"Unable to connect to GLVis server at "
335 << vishost <<
':' << visport << endl;
336 cout <<
"GLVis visualization disabled.\n";
341 sout.precision(precision);
342 sout <<
"solution\n" << *pmesh << u_gf;
347 cout <<
"GLVis visualization paused."
348 <<
" Press space (in the GLVis window) to resume it.\n";
357 cout <<
"Integrating the ODE ..." << endl;
361 ode_solver->
Init(oper);
364 bool last_step =
false;
365 for (
int ti = 1; !last_step; ti++)
367 double dt_real = min(dt, t_final - t);
374 ode_solver->
Step(u, t, dt_real);
376 last_step = (t >= t_final - 1e-8*dt);
378 if (last_step || (ti % vis_steps) == 0)
382 cout <<
"step " << ti <<
", t = " << t << endl;
390 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
391 sout <<
"solution\n" << *pmesh << u_gf << flush;
401 oper.SetParameters(u);
412 ostringstream sol_name;
413 sol_name <<
"ex16-final." << setfill(
'0') << setw(6) << myid;
414 ofstream osol(sol_name.str().c_str());
415 osol.precision(precision);
429 double kap,
const Vector &u)
431 T(NULL), current_dt(0.0),
432 M_solver(f.GetComm()), T_solver(f.GetComm()), z(height)
434 const double rel_tol = 1e-8;
439 M->FormSystemMatrix(ess_tdof_list, Mmat);
441 M_solver.iterative_mode =
false;
442 M_solver.SetRelTol(rel_tol);
443 M_solver.SetAbsTol(0.0);
444 M_solver.SetMaxIter(100);
445 M_solver.SetPrintLevel(0);
446 M_prec.SetType(HypreSmoother::Jacobi);
447 M_solver.SetPreconditioner(M_prec);
448 M_solver.SetOperator(Mmat);
453 T_solver.iterative_mode =
false;
454 T_solver.SetRelTol(rel_tol);
455 T_solver.SetAbsTol(0.0);
456 T_solver.SetMaxIter(100);
457 T_solver.SetPrintLevel(0);
458 T_solver.SetPreconditioner(T_prec);
470 M_solver.Mult(z, du_dt);
473 void ConductionOperator::ImplicitSolve(
const double dt,
481 T =
Add(1.0, Mmat, dt, Kmat);
483 T_solver.SetOperator(*T);
485 MFEM_VERIFY(dt == current_dt,
"");
488 T_solver.Mult(z, du_dt);
491 void ConductionOperator::SundialsSolve(
const double dt,
Vector &b)
494 if (!T || dt != current_dt)
497 T =
Add(1.0, Mmat, dt, Kmat);
499 T_solver.SetOperator(*T);
505 void ConductionOperator::SetParameters(
const Vector &u)
508 u_alpha_gf.SetFromTrueDofs(u);
509 for (
int i = 0; i < u_alpha_gf.Size(); i++)
521 K->FormSystemMatrix(ess_tdof_list, Kmat);
526 ConductionOperator::~ConductionOperator()
534 int SundialsJacSolver::InitSystem(
void *sundials_mem)
539 oper =
dynamic_cast<ConductionOperator*
>(td_oper);
540 MFEM_VERIFY(oper,
"operator is not ConductionOperator");
548 int SundialsJacSolver::SetupSystem(
void *sundials_mem,
int conv_fail,
550 int &jac_cur,
Vector &v_temp1,
558 int SundialsJacSolver::SolveSystem(
void *sundials_mem,
Vector &b,
562 oper->SundialsSolve(GetTimeStep(sundials_mem), b);
567 int SundialsJacSolver::FreeSystem(
void *sundials_mem)
Conjugate gradient method.
double InitialTemperature(const Vector &x)
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Base abstract class for time dependent operators.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
double Norml2() const
Returns the l2 norm of the vector.
virtual void Step(Vector &x, double &t, double &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
void SetStepMode(int itask)
ARKode supports two modes, specified by itask: ARK_NORMAL (default) and ARK_ONE_STEP.
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
virtual void Save()
Save the collection and a VisIt root file.
virtual void Save(std::ostream &out) const
void SetLinearSolver(SundialsODELinearSolver &ls_spec)
Abstract parallel finite element space.
virtual void ProjectCoefficient(Coefficient &coeff)
int main(int argc, char *argv[])
Backward Euler ODE solver. L-stable.
void SetSStolerances(double reltol, double abstol)
Specify the scalar relative and scalar absolute tolerances.
void SetERKTableNum(int table_num)
Choose a specific Butcher table for explicit RK method.
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Wrapper for SUNDIALS' CVODE library – Multi-step time integration.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Data collection with VisIt I/O routines.
HYPRE_Int GlobalTrueVSize() const
void SetMaxStep(double dt_max)
Set the maximum time step of the linear multistep method.
virtual void Print(std::ostream &out=mfem::out) const
Parallel smoothers in hypre.
void PrintUsage(std::ostream &out) const
void SetTime(double t)
Set physical time (for time-dependent simulations)
Wrapper for SUNDIALS' ARKODE library – Runge-Kutta time integration.
void SetLinearSolver(SundialsODELinearSolver &ls_spec)
Set a custom Jacobian system solver for implicit methods.
The classical explicit forth-order Runge-Kutta method, RK4.
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Abstract base class, wrapping the custom linear solvers interface in SUNDIALS' CVODE and ARKODE solve...
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
void PrintOptions(std::ostream &out) const
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
int open(const char hostname[], int port)
class for C-function coefficient
Arbitrary order H1-conforming (continuous) finite elements.
void PrintInfo() const
Print CVODE statistics.
Class for parallel grid function.
The classical forward Euler method.
void SetMaxStep(double dt_max)
Set the maximum time step of the Runge-Kutta method.
Wrapper for hypre's ParCSR matrix class.
Class for parallel meshes.
void PrintInfo() const
Print ARKODE statistics.
void SetStepMode(int itask)
CVode supports two modes, specified by itask: CV_NORMAL (default) and CV_ONE_STEP.