80 virtual void ImplicitSolve(
const double dt,
const Vector &
u,
Vector &k);
95 virtual int SUNImplicitSetup(
const Vector &x,
const Vector &fx,
96 int jok,
int *jcur,
double gamma);
100 virtual int SUNImplicitSolve(
const Vector &
b,
Vector &x,
double tol);
103 void SetParameters(
const Vector &
u);
105 virtual ~ConductionOperator();
110 int main(
int argc,
char *argv[])
116 const char *mesh_file =
"../../data/star.mesh";
119 int ode_solver_type = 9;
120 double t_final = 0.5;
122 double alpha = 1.0e-2;
124 bool visualization =
true;
129 const double reltol = 1e-4, abstol = 1e-4;
132 cout.precision(precision);
135 args.
AddOption(&mesh_file,
"-m",
"--mesh",
136 "Mesh file to use.");
137 args.
AddOption(&ref_levels,
"-r",
"--refine",
138 "Number of times to refine the mesh uniformly.");
140 "Order (degree) of the finite elements.");
141 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
143 "1 - Forward Euler,\n\t" 147 "5 - Backward Euler,\n\t" 150 "8 - CVODE (implicit Adams),\n\t" 151 "9 - CVODE (implicit BDF),\n\t" 152 "10 - ARKODE (default explicit),\n\t" 153 "11 - ARKODE (explicit Fehlberg-6-4-5),\n\t" 154 "12 - ARKODE (default impicit).");
155 args.
AddOption(&t_final,
"-tf",
"--t-final",
156 "Final time; start time is 0.");
157 args.
AddOption(&dt,
"-dt",
"--time-step",
160 "Alpha coefficient.");
162 "Kappa coefficient offset.");
163 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
164 "--no-visualization",
165 "Enable or disable GLVis visualization.");
166 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
167 "--no-visit-datafiles",
168 "Save data files for VisIt (visit.llnl.gov) visualization.");
169 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
170 "Visualize every n-th timestep.");
177 if (ode_solver_type < 1 || ode_solver_type > 12)
179 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
186 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
192 for (
int lev = 0; lev < ref_levels; lev++)
203 cout <<
"Number of temperature unknowns: " << fe_size << endl;
219 ofstream omesh(
"ex16.mesh");
220 omesh.precision(precision);
222 ofstream osol(
"ex16-init.gf");
223 osol.precision(precision);
244 cout <<
"Unable to connect to GLVis server at " 246 visualization =
false;
247 cout <<
"GLVis visualization disabled.\n";
251 sout.precision(precision);
252 sout <<
"solution\n" << *mesh << u_gf;
255 cout <<
"GLVis visualization paused." 256 <<
" Press space (in the GLVis window) to resume it.\n";
265 switch (ode_solver_type)
269 case 2: ode_solver =
new RK2Solver(0.5);
break;
271 case 4: ode_solver =
new RK4Solver;
break;
282 ode_solver = cvode;
break;
288 ode_solver = cvode;
break;
296 if (ode_solver_type == 11)
300 ode_solver = arkode;
break;
306 ode_solver = arkode;
break;
310 if (ode_solver_type < 8) { ode_solver->
Init(oper); }
319 cout <<
"Integrating the ODE ..." << endl;
323 bool last_step =
false;
324 for (
int ti = 1; !last_step; ti++)
326 double dt_real = min(dt, t_final -
t);
333 ode_solver->
Step(
u,
t, dt_real);
335 last_step = (
t >= t_final - 1e-8*dt);
337 if (last_step || (ti % vis_steps) == 0)
339 cout <<
"step " << ti <<
", t = " <<
t << endl;
346 sout <<
"solution\n" << *mesh << u_gf << flush;
356 oper.SetParameters(
u);
364 ofstream osol(
"ex16-final.gf");
365 osol.precision(precision);
381 const double rel_tol = 1e-8;
386 M->FormSystemMatrix(ess_tdof_list, Mmat);
388 M_solver.iterative_mode =
false;
389 M_solver.SetRelTol(rel_tol);
390 M_solver.SetAbsTol(0.0);
391 M_solver.SetMaxIter(50);
392 M_solver.SetPrintLevel(0);
393 M_solver.SetPreconditioner(M_prec);
394 M_solver.SetOperator(Mmat);
399 T_solver.iterative_mode =
false;
400 T_solver.SetRelTol(rel_tol);
401 T_solver.SetAbsTol(0.0);
402 T_solver.SetMaxIter(100);
403 T_solver.SetPrintLevel(0);
404 T_solver.SetPreconditioner(T_prec);
416 M_solver.Mult(z, du_dt);
419 void ConductionOperator::ImplicitSolve(
const double dt,
426 T =
Add(1.0, Mmat, dt, Kmat);
427 T_solver.SetOperator(*T);
430 T_solver.Mult(z, du_dt);
433 void ConductionOperator::SetParameters(
const Vector &
u)
436 u_alpha_gf.SetFromTrueDofs(
u);
437 for (
int i = 0; i < u_alpha_gf.Size(); i++)
449 K->FormSystemMatrix(ess_tdof_list, Kmat);
452 int ConductionOperator::SUNImplicitSetup(
const Vector &x,
453 const Vector &fx,
int jok,
int *jcur,
458 T =
Add(1.0, Mmat, gamma, Kmat);
459 T_solver.SetOperator(*T);
464 int ConductionOperator::SUNImplicitSolve(
const Vector &
b,
Vector &x,
double tol)
472 ConductionOperator::~ConductionOperator()
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
Conjugate gradient method.
Class for grid function - Vector with associated FE space.
double InitialTemperature(const Vector &x)
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
Data type for scaled Jacobi-type smoother of sparse matrix.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
void PrintOptions(std::ostream &out) const
Print the options.
int Dimension() const
Dimension of the reference space used within the elements.
void SetStepMode(int itask)
Select the ARKode step mode: ARK_NORMAL (default) or ARK_ONE_STEP.
void SetERKTableNum(ARKODE_ERKTableID table_id)
Choose a specific Butcher table for an explicit RK method.
void PrintUsage(std::ostream &out) const
Print the usage message.
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)
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]...
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
void SetMaxStep(double dt_max)
Set the maximum time step.
double RealTime()
Return the number of real seconds elapsed since the stopwatch was started.
bool Good() const
Return true if the command line options were parsed successfully.
Backward Euler ODE solver. L-stable.
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.
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
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.
void PrintInfo() const
Print various ARKStep statistics.
void SetTime(double t)
Set physical time (for time-dependent simulations)
void Start()
Start the stopwatch. The elapsed time is not cleared.
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
The classical explicit forth-order Runge-Kutta method, RK4.
void PrintInfo() const
Print various CVODE statistics.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual void Save() override
Save the collection and a VisIt root file.
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.
constexpr ARKODE_ERKTableID ARKODE_FEHLBERG_13_7_8
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
double Norml2() const
Returns the l2 norm of the vector.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
A general function coefficient.
virtual void Print(std::ostream &os=mfem::out) const
Arbitrary order H1-conforming (continuous) finite elements.
double u(const Vector &xvec)
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
The classical forward Euler method.
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults...
int main(int argc, char *argv[])
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)
virtual void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
void SetStepMode(int itask)
Select the CVODE step mode: CV_NORMAL (default) or CV_ONE_STEP.