85 void SetParameters(
const Vector &
u);
87 ~ConductionOperator()
override;
92int main(
int argc,
char *argv[])
95 const char *mesh_file =
"../data/star.mesh";
99 int ode_solver_type = 23;
105 bool visualization =
true;
110 cout.precision(precision);
113 args.
AddOption(&mesh_file,
"-m",
"--mesh",
114 "Mesh file to use.");
115 args.
AddOption(&ref_levels,
"-r",
"--refine",
116 "Number of times to refine the mesh uniformly.");
118 "Order (degree) of the finite elements.");
119 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
121 args.
AddOption(&t_final,
"-tf",
"--t-final",
122 "Final time; start time is 0.");
123 args.
AddOption(&dt,
"-dt",
"--time-step",
126 "Alpha coefficient.");
128 "Kappa coefficient offset.");
129 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
130 "--no-visualization",
131 "Enable or disable GLVis visualization.");
132 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
133 "--no-visit-datafiles",
134 "Save data files for VisIt (visit.llnl.gov) visualization.");
135 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
136 "Visualize every n-th timestep.");
147 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
158 for (
int lev = 0; lev < ref_levels; lev++)
169 cout <<
"Number of temperature unknowns: " << fe_size << endl;
185 ofstream omesh(
"ex16.mesh");
186 omesh.precision(precision);
188 ofstream osol(
"ex16-init.gf");
189 osol.precision(precision);
210 cout <<
"Unable to connect to GLVis server at "
211 <<
vishost <<
':' << visport << endl;
212 visualization =
false;
213 cout <<
"GLVis visualization disabled.\n";
217 sout.precision(precision);
218 sout <<
"solution\n" << *mesh << u_gf;
221 cout <<
"GLVis visualization paused."
222 <<
" Press space (in the GLVis window) to resume it.\n";
228 ode_solver->Init(oper);
231 bool last_step =
false;
232 for (
int ti = 1; !last_step; ti++)
234 if (t + dt >= t_final - dt/2)
239 ode_solver->Step(
u, t, dt);
241 if (last_step || (ti % vis_steps) == 0)
243 cout <<
"step " << ti <<
", t = " << t << endl;
248 sout <<
"solution\n" << *mesh << u_gf << flush;
258 oper.SetParameters(
u);
264 ofstream osol(
"ex16-final.gf");
265 osol.precision(precision);
278 M(NULL), K(NULL), T(NULL), current_dt(0.0), z(height)
280 const real_t rel_tol = 1e-8;
285 M->FormSystemMatrix(ess_tdof_list, Mmat);
287 M_solver.iterative_mode =
false;
288 M_solver.SetRelTol(rel_tol);
289 M_solver.SetAbsTol(0.0);
290 M_solver.SetMaxIter(30);
291 M_solver.SetPrintLevel(0);
292 M_solver.SetPreconditioner(M_prec);
293 M_solver.SetOperator(Mmat);
298 T_solver.iterative_mode =
false;
299 T_solver.SetRelTol(rel_tol);
300 T_solver.SetAbsTol(0.0);
301 T_solver.SetMaxIter(100);
302 T_solver.SetPrintLevel(0);
303 T_solver.SetPreconditioner(T_prec);
308void ConductionOperator::Mult(
const Vector &
u,
Vector &du_dt)
const
315 M_solver.
Mult(z, du_dt);
318void ConductionOperator::ImplicitSolve(
const real_t dt,
326 T =
Add(1.0, Mmat, dt, Kmat);
330 MFEM_VERIFY(dt == current_dt,
"");
333 T_solver.
Mult(z, du_dt);
336void ConductionOperator::SetParameters(
const Vector &
u)
339 u_alpha_gf.SetFromTrueDofs(
u);
340 for (
int i = 0; i < u_alpha_gf.Size(); i++)
342 u_alpha_gf(i) = kappa + alpha*u_alpha_gf(i);
357ConductionOperator::~ConductionOperator()
Conjugate gradient method.
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
void SetOperator(const Operator &op) override
Set/update the solver for the given operator.
Data type for scaled Jacobi-type smoother of sparse matrix.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetTime(real_t t)
Set physical time (for time-dependent simulations)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
A general function coefficient.
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
Class for grid function - Vector with associated FE space.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
Arbitrary order H1-conforming (continuous) finite elements.
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
int Dimension() const
Dimension of the reference space used within the elements.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
static MFEM_EXPORT std::string Types
static MFEM_EXPORT std::unique_ptr< ODESolver > Select(const int ode_solver_type)
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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...
bool Good() const
Return true if the command line options were parsed successfully.
void Mult(const Vector &x, Vector &y) const override
Matrix vector multiplication.
Base abstract class for first order time dependent operators.
void Neg()
(*this) = -(*this)
real_t Norml2() const
Returns the l2 norm of the vector.
Data collection with VisIt I/O routines.
void Save() override
Save the collection and a VisIt root file.
void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
real_t InitialTemperature(const Vector &x)
real_t u(const Vector &xvec)
std::function< real_t(const Vector &)> f(real_t mass_coeff)
void Add(const DenseMatrix &A, const DenseMatrix &B, real_t alpha, DenseMatrix &C)
C = A + alpha*B.