79 virtual void ImplicitSolve(
const double dt,
const Vector &u,
Vector &k);
82 void SetParameters(
const Vector &u);
84 virtual ~ConductionOperator();
89 int main(
int argc,
char *argv[])
93 MPI_Init(&argc, &argv);
94 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
95 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
98 const char *mesh_file =
"../data/star.mesh";
99 int ser_ref_levels = 2;
100 int par_ref_levels = 1;
102 int ode_solver_type = 3;
103 double t_final = 0.5;
105 double alpha = 1.0e-2;
107 bool visualization =
true;
112 cout.precision(precision);
115 args.
AddOption(&mesh_file,
"-m",
"--mesh",
116 "Mesh file to use.");
117 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
118 "Number of times to refine the mesh uniformly in serial.");
119 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
120 "Number of times to refine the mesh uniformly in parallel.");
122 "Order (degree) of the finite elements.");
123 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
124 "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t"
125 "\t 11 - Forward Euler, 12 - RK2, 13 - RK3 SSP, 14 - RK4.");
126 args.
AddOption(&t_final,
"-tf",
"--t-final",
127 "Final time; start time is 0.");
128 args.
AddOption(&dt,
"-dt",
"--time-step",
131 "Alpha coefficient.");
133 "Kappa coefficient offset.");
134 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
135 "--no-visualization",
136 "Enable or disable GLVis visualization.");
137 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
138 "--no-visit-datafiles",
139 "Save data files for VisIt (visit.llnl.gov) visualization.");
140 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
141 "Visualize every n-th timestep.");
158 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
165 switch (ode_solver_type)
173 case 12: ode_solver =
new RK2Solver(0.5);
break;
175 case 14: ode_solver =
new RK4Solver;
break;
181 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
188 for (
int lev = 0; lev < ser_ref_levels; lev++)
198 for (
int lev = 0; lev < par_ref_levels; lev++)
211 cout <<
"Number of temperature unknowns: " << fe_size << endl;
224 ConductionOperator oper(fespace, alpha, kappa, u);
228 ostringstream mesh_name, sol_name;
229 mesh_name <<
"ex16-mesh." << setfill(
'0') << setw(6) << myid;
230 sol_name <<
"ex16-init." << setfill(
'0') << setw(6) << myid;
231 ofstream omesh(mesh_name.str().c_str());
232 omesh.precision(precision);
234 ofstream osol(sol_name.str().c_str());
235 osol.precision(precision);
251 char vishost[] =
"localhost";
253 sout.
open(vishost, visport);
254 sout <<
"parallel " << num_procs <<
" " << myid << endl;
255 int good = sout.good(), all_good;
256 MPI_Allreduce(&good, &all_good, 1, MPI_INT, MPI_MIN, pmesh->
GetComm());
260 visualization =
false;
263 cout <<
"Unable to connect to GLVis server at "
264 << vishost <<
':' << visport << endl;
265 cout <<
"GLVis visualization disabled.\n";
270 sout.precision(precision);
271 sout <<
"solution\n" << *pmesh << u_gf;
276 cout <<
"GLVis visualization paused."
277 <<
" Press space (in the GLVis window) to resume it.\n";
284 ode_solver->
Init(oper);
287 bool last_step =
false;
288 for (
int ti = 1; !last_step; ti++)
290 if (t + dt >= t_final - dt/2)
295 ode_solver->
Step(u, t, dt);
297 if (last_step || (ti % vis_steps) == 0)
301 cout <<
"step " << ti <<
", t = " << t << endl;
307 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
308 sout <<
"solution\n" << *pmesh << u_gf << flush;
318 oper.SetParameters(u);
324 ostringstream sol_name;
325 sol_name <<
"ex16-final." << setfill(
'0') << setw(6) << myid;
326 ofstream osol(sol_name.str().c_str());
327 osol.precision(precision);
341 double kap,
const Vector &u)
343 T(NULL), current_dt(0.0),
344 M_solver(f.GetComm()), T_solver(f.GetComm()), z(height)
346 const double rel_tol = 1e-8;
351 M->FormSystemMatrix(ess_tdof_list, Mmat);
353 M_solver.iterative_mode =
false;
354 M_solver.SetRelTol(rel_tol);
355 M_solver.SetAbsTol(0.0);
356 M_solver.SetMaxIter(100);
357 M_solver.SetPrintLevel(0);
358 M_prec.SetType(HypreSmoother::Jacobi);
359 M_solver.SetPreconditioner(M_prec);
360 M_solver.SetOperator(Mmat);
365 T_solver.iterative_mode =
false;
366 T_solver.SetRelTol(rel_tol);
367 T_solver.SetAbsTol(0.0);
368 T_solver.SetMaxIter(100);
369 T_solver.SetPrintLevel(0);
370 T_solver.SetPreconditioner(T_prec);
382 M_solver.Mult(z, du_dt);
385 void ConductionOperator::ImplicitSolve(
const double dt,
393 T =
Add(1.0, Mmat, dt, Kmat);
395 T_solver.SetOperator(*T);
397 MFEM_VERIFY(dt == current_dt,
"");
400 T_solver.Mult(z, du_dt);
403 void ConductionOperator::SetParameters(
const Vector &u)
406 u_alpha_gf.SetFromTrueDofs(u);
407 for (
int i = 0; i < u_alpha_gf.Size(); i++)
419 K->FormSystemMatrix(ess_tdof_list, Kmat);
424 ConductionOperator::~ConductionOperator()
Conjugate gradient method.
double InitialTemperature(const Vector &x)
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
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]...
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
Abstract parallel finite element space.
virtual void ProjectCoefficient(Coefficient &coeff)
Backward Euler ODE solver. L-stable.
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 UniformRefinement(int i, const DSTable &, int *, int *, int *)
Data collection with VisIt I/O routines.
HYPRE_Int GlobalTrueVSize() const
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)
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.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
Implicit midpoint method. A-stable, not L-stable.
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.
Class for parallel grid function.
The classical forward Euler method.
Wrapper for hypre's ParCSR matrix class.
Class for parallel meshes.