84 virtual void ImplicitSolve(
const double dt,
const Vector &
u,
Vector &k);
87 void SetParameters(
const Vector &
u);
89 virtual ~ConductionOperator();
94 int main(
int argc,
char *argv[])
97 Mpi::Init(argc, argv);
98 int num_procs = Mpi::WorldSize();
99 int myid = Mpi::WorldRank();
103 const char *mesh_file =
"../data/star.mesh";
104 int ser_ref_levels = 2;
105 int par_ref_levels = 1;
107 int ode_solver_type = 3;
108 double t_final = 0.5;
110 double alpha = 1.0e-2;
112 bool visualization =
true;
118 cout.precision(precision);
121 args.
AddOption(&mesh_file,
"-m",
"--mesh",
122 "Mesh file to use.");
123 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
124 "Number of times to refine the mesh uniformly in serial.");
125 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
126 "Number of times to refine the mesh uniformly in parallel.");
128 "Order (degree) of the finite elements.");
129 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
130 "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t" 131 "\t 11 - Forward Euler, 12 - RK2, 13 - RK3 SSP, 14 - RK4.");
132 args.
AddOption(&t_final,
"-tf",
"--t-final",
133 "Final time; start time is 0.");
134 args.
AddOption(&dt,
"-dt",
"--time-step",
137 "Alpha coefficient.");
139 "Kappa coefficient offset.");
140 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
141 "--no-visualization",
142 "Enable or disable GLVis visualization.");
143 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
144 "--no-visit-datafiles",
145 "Save data files for VisIt (visit.llnl.gov) visualization.");
146 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
147 "Visualize every n-th timestep.");
148 args.
AddOption(&adios2,
"-adios2",
"--adios2-streams",
"-no-adios2",
149 "--no-adios2-streams",
150 "Save data using adios2 streams.");
166 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
173 switch (ode_solver_type)
181 case 12: ode_solver =
new RK2Solver(0.5);
break;
183 case 14: ode_solver =
new RK4Solver;
break;
190 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
198 for (
int lev = 0; lev < ser_ref_levels; lev++)
208 for (
int lev = 0; lev < par_ref_levels; lev++)
221 cout <<
"Number of temperature unknowns: " << fe_size << endl;
238 ostringstream mesh_name, sol_name;
239 mesh_name <<
"ex16-mesh." << setfill(
'0') << setw(6) << myid;
240 sol_name <<
"ex16-init." << setfill(
'0') << setw(6) << myid;
241 ofstream omesh(mesh_name.str().c_str());
242 omesh.precision(precision);
244 ofstream osol(sol_name.str().c_str());
245 osol.precision(precision);
260 #ifdef MFEM_USE_ADIOS2 264 std::string postfix(mesh_file);
265 postfix.erase(0, std::string(
"../data/").size() );
266 postfix +=
"_o" + std::to_string(order);
267 postfix +=
"_solver" + std::to_string(ode_solver_type);
268 const std::string collection_name =
"ex16-p-" + postfix +
".bp";
271 adios2_dc->
SetParameter(
"SubStreams", std::to_string(num_procs/2) );
285 sout <<
"parallel " << num_procs <<
" " << myid << endl;
286 int good = sout.good(), all_good;
287 MPI_Allreduce(&good, &all_good, 1, MPI_INT, MPI_MIN, pmesh->
GetComm());
291 visualization =
false;
294 cout <<
"Unable to connect to GLVis server at " 296 cout <<
"GLVis visualization disabled.\n";
301 sout.precision(precision);
302 sout <<
"solution\n" << *pmesh << u_gf;
307 cout <<
"GLVis visualization paused." 308 <<
" Press space (in the GLVis window) to resume it.\n";
315 ode_solver->
Init(oper);
318 bool last_step =
false;
319 for (
int ti = 1; !last_step; ti++)
321 if (
t + dt >= t_final - dt/2)
326 ode_solver->
Step(
u,
t, dt);
328 if (last_step || (ti % vis_steps) == 0)
332 cout <<
"step " << ti <<
", t = " <<
t << endl;
338 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
339 sout <<
"solution\n" << *pmesh << u_gf << flush;
349 #ifdef MFEM_USE_ADIOS2 358 oper.SetParameters(
u);
361 #ifdef MFEM_USE_ADIOS2 371 ostringstream sol_name;
372 sol_name <<
"ex16-final." << setfill(
'0') << setw(6) << myid;
373 ofstream osol(sol_name.str().c_str());
388 T(NULL), current_dt(0.0),
389 M_solver(
f.GetComm()), T_solver(
f.GetComm()), z(height)
391 const double rel_tol = 1e-8;
396 M->FormSystemMatrix(ess_tdof_list, Mmat);
398 M_solver.iterative_mode =
false;
399 M_solver.SetRelTol(rel_tol);
400 M_solver.SetAbsTol(0.0);
401 M_solver.SetMaxIter(100);
402 M_solver.SetPrintLevel(0);
403 M_prec.SetType(HypreSmoother::Jacobi);
404 M_solver.SetPreconditioner(M_prec);
405 M_solver.SetOperator(Mmat);
410 T_solver.iterative_mode =
false;
411 T_solver.SetRelTol(rel_tol);
412 T_solver.SetAbsTol(0.0);
413 T_solver.SetMaxIter(100);
414 T_solver.SetPrintLevel(0);
415 T_solver.SetPreconditioner(T_prec);
427 M_solver.Mult(z, du_dt);
430 void ConductionOperator::ImplicitSolve(
const double dt,
438 T =
Add(1.0, Mmat, dt, Kmat);
440 T_solver.SetOperator(*T);
442 MFEM_VERIFY(dt == current_dt,
"");
445 T_solver.Mult(z, du_dt);
448 void ConductionOperator::SetParameters(
const Vector &
u)
451 u_alpha_gf.SetFromTrueDofs(
u);
452 for (
int i = 0; i < u_alpha_gf.Size(); i++)
464 K->FormSystemMatrix(ess_tdof_list, Kmat);
469 ConductionOperator::~ConductionOperator()
Conjugate gradient method.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void PrintOptions(std::ostream &out) const
Print the options.
void PrintUsage(std::ostream &out) const
Print the usage message.
Base abstract class for first order time dependent operators.
int main(int argc, char *argv[])
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)
bool Good() const
Return true if the command line options were parsed successfully.
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 Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
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 ...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Data collection with VisIt I/O routines.
void SetParameter(const std::string key, const std::string value) noexcept
Parallel smoothers in hypre.
void SetTime(double t)
Set physical time (for time-dependent simulations)
HYPRE_BigInt GlobalTrueVSize() const
The classical explicit forth-order Runge-Kutta method, RK4.
int precision
Precision (number of digits) used for the text output of doubles.
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.
virtual void Save(std::ostream &out) const
Implicit midpoint method. A-stable, not L-stable.
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.
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Arbitrary order H1-conforming (continuous) finite elements.
void Print(std::ostream &out=mfem::out) const override
double u(const Vector &xvec)
Class for parallel grid function.
double InitialTemperature(const Vector &x)
The classical forward Euler method.
Wrapper for hypre's ParCSR matrix class.
Class for parallel meshes.
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
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.