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[])
106 MPI_Init(&argc, &argv);
107 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
108 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
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.");
183 if (ode_solver_type < 1 || ode_solver_type > 12)
187 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
196 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
202 for (
int lev = 0; lev < ser_ref_levels; lev++)
212 for (
int lev = 0; lev < par_ref_levels; lev++)
225 cout <<
"Number of temperature unknowns: " << fe_size << endl;
238 ConductionOperator oper(fespace, alpha, kappa, u);
242 ostringstream mesh_name, sol_name;
243 mesh_name <<
"ex16-mesh." << setfill(
'0') << setw(6) << myid;
244 sol_name <<
"ex16-init." << setfill(
'0') << setw(6) << myid;
245 ofstream omesh(mesh_name.str().c_str());
246 omesh.precision(precision);
248 ofstream osol(sol_name.str().c_str());
249 osol.precision(precision);
267 sout.
open(vishost, visport);
268 sout <<
"parallel " << num_procs <<
" " << myid << endl;
269 int good = sout.good(), all_good;
270 MPI_Allreduce(&good, &all_good, 1, MPI_INT, MPI_MIN, pmesh->
GetComm());
274 visualization =
false;
277 cout <<
"Unable to connect to GLVis server at "
278 << vishost <<
':' << visport << endl;
279 cout <<
"GLVis visualization disabled.\n";
284 sout.precision(precision);
285 sout <<
"solution\n" << *pmesh << u_gf;
290 cout <<
"GLVis visualization paused."
291 <<
" Press space (in the GLVis window) to resume it.\n";
301 switch (ode_solver_type)
305 case 2: ode_solver =
new RK2Solver(0.5);
break;
307 case 4: ode_solver =
new RK4Solver;
break;
318 ode_solver = cvode;
break;
324 ode_solver = cvode;
break;
328 arkode =
new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::EXPLICIT);
332 if (ode_solver_type == 11) { arkode->
SetERKTableNum(FEHLBERG_13_7_8); }
333 ode_solver = arkode;
break;
335 arkode =
new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::IMPLICIT);
339 ode_solver = arkode;
break;
343 if (ode_solver_type < 8) { ode_solver->
Init(oper); }
354 cout <<
"Integrating the ODE ..." << endl;
359 bool last_step =
false;
360 for (
int ti = 1; !last_step; ti++)
362 double dt_real = min(dt, t_final - t);
369 ode_solver->
Step(u, t, dt_real);
371 last_step = (t >= t_final - 1e-8*dt);
373 if (last_step || (ti % vis_steps) == 0)
377 cout <<
"step " << ti <<
", t = " << t << endl;
385 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
386 sout <<
"solution\n" << *pmesh << u_gf << flush;
396 oper.SetParameters(u);
407 ostringstream sol_name;
408 sol_name <<
"ex16-final." << setfill(
'0') << setw(6) << myid;
409 ofstream osol(sol_name.str().c_str());
410 osol.precision(precision);
424 double kap,
const Vector &u)
427 M_solver(f.GetComm()), T_solver(f.GetComm()), z(height)
429 const double rel_tol = 1e-8;
434 M->FormSystemMatrix(ess_tdof_list, Mmat);
436 M_solver.iterative_mode =
false;
437 M_solver.SetRelTol(rel_tol);
438 M_solver.SetAbsTol(0.0);
439 M_solver.SetMaxIter(100);
440 M_solver.SetPrintLevel(0);
441 M_prec.SetType(HypreSmoother::Jacobi);
442 M_solver.SetPreconditioner(M_prec);
443 M_solver.SetOperator(Mmat);
448 T_solver.iterative_mode =
false;
449 T_solver.SetRelTol(rel_tol);
450 T_solver.SetAbsTol(0.0);
451 T_solver.SetMaxIter(100);
452 T_solver.SetPrintLevel(0);
453 T_solver.SetPreconditioner(T_prec);
465 M_solver.Mult(z, du_dt);
468 void ConductionOperator::ImplicitSolve(
const double dt,
475 T =
Add(1.0, Mmat, dt, Kmat);
476 T_solver.SetOperator(*T);
479 T_solver.Mult(z, du_dt);
482 int ConductionOperator::SUNImplicitSetup(
const Vector &x,
483 const Vector &fx,
int jok,
int *jcur,
488 T =
Add(1.0, Mmat, gamma, Kmat);
489 T_solver.SetOperator(*T);
494 int ConductionOperator::SUNImplicitSolve(
const Vector &
b,
Vector &x,
double tol)
502 void ConductionOperator::SetParameters(
const Vector &u)
505 u_alpha_gf.SetFromTrueDofs(u);
506 for (
int i = 0; i < u_alpha_gf.Size(); i++)
518 K->FormSystemMatrix(ess_tdof_list, Kmat);
521 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]...
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
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)
int main(int argc, char *argv[])
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.
HYPRE_Int GlobalTrueVSize() const
void SetMaxStep(double dt_max)
Set the maximum time step.
virtual void Print(std::ostream &out=mfem::out) const
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()
Clear the elapsed time and start the stopwatch.
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.
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.
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.