83 virtual void ImplicitSolve(
const double dt,
const Vector &
u,
Vector &k);
87 virtual int SUNImplicitSetup(
const Vector &x,
const Vector &fx,
88 int jok,
int *jcur,
double gamma);
92 virtual int SUNImplicitSolve(
const Vector &
b,
Vector &x,
double tol);
95 void SetParameters(
const Vector &
u);
97 virtual ~ConductionOperator();
102 int main(
int argc,
char *argv[])
105 Mpi::Init(argc, argv);
106 int num_procs = Mpi::WorldSize();
107 int myid = Mpi::WorldRank();
111 const char *mesh_file =
"../../data/star.mesh";
112 int ser_ref_levels = 2;
113 int par_ref_levels = 1;
115 int ode_solver_type = 9;
116 double t_final = 0.5;
118 double alpha = 1.0e-2;
120 bool visualization =
true;
125 const double reltol = 1e-4, abstol = 1e-4;
128 cout.precision(precision);
131 args.
AddOption(&mesh_file,
"-m",
"--mesh",
132 "Mesh file to use.");
133 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
134 "Number of times to refine the mesh uniformly in serial.");
135 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
136 "Number of times to refine the mesh uniformly in parallel.");
138 "Order (degree) of the finite elements.");
139 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
141 "1 - Forward Euler,\n\t"
145 "5 - Backward Euler,\n\t"
148 "8 - CVODE (implicit Adams),\n\t"
149 "9 - CVODE (implicit BDF),\n\t"
150 "10 - ARKODE (default explicit),\n\t"
151 "11 - ARKODE (explicit Fehlberg-6-4-5),\n\t"
152 "12 - ARKODE (default impicit).");
153 args.
AddOption(&t_final,
"-tf",
"--t-final",
154 "Final time; start time is 0.");
155 args.
AddOption(&dt,
"-dt",
"--time-step",
158 "Alpha coefficient.");
160 "Kappa coefficient offset.");
161 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
162 "--no-visualization",
163 "Enable or disable GLVis visualization.");
164 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
165 "--no-visit-datafiles",
166 "Save data files for VisIt (visit.llnl.gov) visualization.");
167 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
168 "Visualize every n-th timestep.");
182 if (ode_solver_type < 1 || ode_solver_type > 12)
186 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
194 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
200 for (
int lev = 0; lev < ser_ref_levels; lev++)
210 for (
int lev = 0; lev < par_ref_levels; lev++)
223 cout <<
"Number of temperature unknowns: " << fe_size << endl;
236 ConductionOperator oper(fespace, alpha, kappa, u);
240 ostringstream mesh_name, sol_name;
241 mesh_name <<
"ex16-mesh." << setfill(
'0') << setw(6) << myid;
242 sol_name <<
"ex16-init." << setfill(
'0') << setw(6) << myid;
243 ofstream omesh(mesh_name.str().c_str());
244 omesh.precision(precision);
246 ofstream osol(sol_name.str().c_str());
247 osol.precision(precision);
265 sout.
open(vishost, visport);
266 sout <<
"parallel " << num_procs <<
" " << myid << endl;
267 int good = sout.good(), all_good;
268 MPI_Allreduce(&good, &all_good, 1, MPI_INT, MPI_MIN, pmesh->
GetComm());
272 visualization =
false;
275 cout <<
"Unable to connect to GLVis server at "
276 << vishost <<
':' << visport << endl;
277 cout <<
"GLVis visualization disabled.\n";
282 sout.precision(precision);
283 sout <<
"solution\n" << *pmesh << u_gf;
288 cout <<
"GLVis visualization paused."
289 <<
" Press space (in the GLVis window) to resume it.\n";
299 switch (ode_solver_type)
303 case 2: ode_solver =
new RK2Solver(0.5);
break;
305 case 4: ode_solver =
new RK4Solver;
break;
316 ode_solver = cvode;
break;
322 ode_solver = cvode;
break;
326 arkode =
new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::EXPLICIT);
330 if (ode_solver_type == 11) { arkode->
SetERKTableNum(FEHLBERG_13_7_8); }
331 ode_solver = arkode;
break;
333 arkode =
new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::IMPLICIT);
337 ode_solver = arkode;
break;
341 if (ode_solver_type < 8) { ode_solver->
Init(oper); }
352 cout <<
"Integrating the ODE ..." << endl;
357 bool last_step =
false;
358 for (
int ti = 1; !last_step; ti++)
360 double dt_real = min(dt, t_final - t);
367 ode_solver->
Step(u, t, dt_real);
369 last_step = (t >= t_final - 1e-8*dt);
371 if (last_step || (ti % vis_steps) == 0)
375 cout <<
"step " << ti <<
", t = " << t << endl;
383 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
384 sout <<
"solution\n" << *pmesh << u_gf << flush;
394 oper.SetParameters(u);
405 ostringstream sol_name;
406 sol_name <<
"ex16-final." << setfill(
'0') << setw(6) << myid;
407 ofstream osol(sol_name.str().c_str());
408 osol.precision(precision);
423 M_solver(f.GetComm()), T_solver(f.GetComm()), z(height)
425 const double rel_tol = 1e-8;
430 M->FormSystemMatrix(ess_tdof_list, Mmat);
432 M_solver.iterative_mode =
false;
433 M_solver.SetRelTol(rel_tol);
434 M_solver.SetAbsTol(0.0);
435 M_solver.SetMaxIter(100);
436 M_solver.SetPrintLevel(0);
437 M_prec.SetType(HypreSmoother::Jacobi);
438 M_solver.SetPreconditioner(M_prec);
439 M_solver.SetOperator(Mmat);
444 T_solver.iterative_mode =
false;
445 T_solver.SetRelTol(rel_tol);
446 T_solver.SetAbsTol(0.0);
447 T_solver.SetMaxIter(100);
448 T_solver.SetPrintLevel(0);
449 T_solver.SetPreconditioner(T_prec);
461 M_solver.Mult(z, du_dt);
464 void ConductionOperator::ImplicitSolve(
const double dt,
471 T =
Add(1.0, Mmat, dt, Kmat);
472 T_solver.SetOperator(*T);
475 T_solver.Mult(z, du_dt);
478 int ConductionOperator::SUNImplicitSetup(
const Vector &x,
479 const Vector &fx,
int jok,
int *jcur,
484 T =
Add(1.0, Mmat, gamma, Kmat);
485 T_solver.SetOperator(*T);
490 int ConductionOperator::SUNImplicitSolve(
const Vector &
b,
Vector &x,
double tol)
498 void ConductionOperator::SetParameters(
const Vector &u)
501 u_alpha_gf.SetFromTrueDofs(u);
502 for (
int i = 0; i < u_alpha_gf.Size(); i++)
514 K->FormSystemMatrix(ess_tdof_list, Kmat);
517 ConductionOperator::~ConductionOperator()
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
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.
void SetStepMode(int itask)
Select the ARKode step mode: ARK_NORMAL (default) or ARK_ONE_STEP.
Base abstract class for first order 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]...
HYPRE_BigInt GlobalTrueVSize() const
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 SetMaxStep(double dt_max)
Set the maximum time step.
Abstract parallel finite element space.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Backward Euler ODE solver. L-stable.
int close()
Close the socketstream.
void Stop()
Stop the stopwatch.
Interface to ARKode's ARKStep module – additive Runge-Kutta methods.
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.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Interface to the CVODE library – linear multi-step methods.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Data collection with VisIt I/O routines.
void SetMaxStep(double dt_max)
Set the maximum time step.
Parallel smoothers in hypre.
void PrintUsage(std::ostream &out) const
Print the usage message.
void SetTime(double t)
Set physical time (for time-dependent simulations)
void Start()
Start the stopwatch. The elapsed time is not cleared.
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)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
Third-order, strong stability preserving (SSP) Runge-Kutta method.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
void PrintInfo() const
Print various ARKStep statistics.
void PrintOptions(std::ostream &out) const
Print the options.
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
A general function coefficient.
Arbitrary order H1-conforming (continuous) finite elements.
void PrintInfo() const
Print various CVODE statistics.
void Print(std::ostream &out=mfem::out) const override
double u(const Vector &xvec)
Class for parallel grid function.
void SetERKTableNum(int table_num)
Choose a specific Butcher table for an explicit RK method.
The classical forward Euler method.
Wrapper for hypre's ParCSR matrix class.
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults...
Class for parallel meshes.
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.
void Clear()
Clear the elapsed time on the stopwatch and restart it if it's running.
double f(const Vector &p)
void SetStepMode(int itask)
Select the CVODE step mode: CV_NORMAL (default) or CV_ONE_STEP.
bool Good() const
Return true if the command line options were parsed successfully.