libMesh
claw_system.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 #ifndef CLAW_SYSTEM_H
19 #define CLAW_SYSTEM_H
20 
21 // libMesh includes
22 #include "libmesh/linear_implicit_system.h" // base class
23 
24 namespace libMesh
25 {
26 
37 {
38 public:
39 
46 
52  const std::string & name,
53  const unsigned int number);
54 
58  virtual ~ClawSystem ();
59 
64 
68  virtual std::string system_type () const override;
69 
73  static TemporalDiscretizationType string_to_enum(std::string string_type);
74  static std::string enum_to_string(TemporalDiscretizationType enum_type);
75 
80  virtual void process_parameters_file (const std::string& parameters_filename);
81 
86 
91  virtual void assemble_claw_rhs(NumericVector<Number> &) = 0;
92 
96  void assemble_all_matrices();
97 
105  SparseMatrix<Number> & get_avg_matrix(unsigned int dim);
108 
112  void set_delta_t(Real delta_t_in);
113  Real get_delta_t();
114 
118  void set_n_time_steps(unsigned int n_time_steps_in);
119  unsigned int get_n_time_steps();
120 
125  void set_LxF_constant(Real LxF_constant_in);
127 
133 
137  void set_write_interval(unsigned int write_interval_in);
138  unsigned int get_write_interval();
139 
143  void set_time(Real time_in);
144  Real get_time();
145 
149  virtual void print_info();
150 
155 
159  virtual void init_data() override;
160 
161 private:
162 
166  void assemble_mass_matrix();
171 
175  std::unique_ptr<SparseMatrix<Number>> _mass_matrix;
176 
180  std::vector<std::unique_ptr<SparseMatrix<Number>>> _advection_matrices;
181 
185  std::vector<std::unique_ptr<SparseMatrix<Number>>> _avg_matrices;
186 
190  std::unique_ptr<SparseMatrix<Number>> _jump_matrix;
191 
196  std::vector<std::unique_ptr<SparseMatrix<Number>>> _boundary_condition_matrices;
197 
203 
208 
212  unsigned int _n_time_steps;
213 
218 
222  unsigned int _write_interval;
223 
228 };
229 
230 } // namespace libMesh
231 
232 #endif
std::vector< std::unique_ptr< SparseMatrix< Number > > > _avg_matrices
The "flux average" element boundary matrices.
Definition: claw_system.h:185
void set_n_time_steps(unsigned int n_time_steps_in)
Set/get the number of time-steps.
Definition: claw_system.C:870
std::unique_ptr< SparseMatrix< Number > > _jump_matrix
The "jump" element boundary matrix.
Definition: claw_system.h:190
This is the EquationSystems class.
SparseMatrix< Number > & get_boundary_condition_matrix(unsigned int dim)
Definition: claw_system.C:848
void assemble_jump_coupling_matrix()
Definition: claw_system.C:590
This class encapsulates functionality that allows us to solve conservation laws.
Definition: claw_system.h:36
TemporalDiscretizationType get_temporal_discretization_type()
Definition: claw_system.C:900
unsigned int dim
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
Real _time
The current time in the conservation law solve.
Definition: claw_system.h:227
LinearImplicitSystem Parent
The type of the parent.
Definition: claw_system.h:63
void assemble_all_matrices()
Assemble the matrices we need to solve a conservation law.
Definition: claw_system.C:262
void set_time(Real time_in)
Set/get the current time in the conservation law solve.
Definition: claw_system.C:918
void assemble_mass_matrix()
Helper functions called by assemble_all_matrices().
Definition: claw_system.C:271
The libMesh namespace provides an interface to certain functionality in the library.
virtual void init_data() override
Initialize the system (e.g.
Definition: claw_system.C:944
void set_temporal_discretization_type(TemporalDiscretizationType td_in)
Set/get the temporal discretization type.
Definition: claw_system.C:894
void set_write_interval(unsigned int write_interval_in)
Set/get write_interval.
Definition: claw_system.C:906
SparseMatrix< Number > & get_advection_matrix(unsigned int dim)
Definition: claw_system.C:822
TemporalDiscretizationType
Define an enumeration for temporal discretization types ForwardEuler = 1st-order Explicit Euler RK4 =...
Definition: claw_system.h:45
unsigned int number() const
Definition: system.h:2350
virtual void process_parameters_file(const std::string &parameters_filename)
Set parameters for this system (e.g.
Definition: claw_system.C:231
Real _delta_t
The time step size.
Definition: claw_system.h:207
SparseMatrix< Number > & get_avg_matrix(unsigned int dim)
Definition: claw_system.C:832
void assemble_boundary_condition_matrices()
Definition: claw_system.C:736
ClawSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
Definition: claw_system.C:39
static std::string enum_to_string(TemporalDiscretizationType enum_type)
Definition: claw_system.C:80
unsigned int get_n_time_steps()
Definition: claw_system.C:876
SparseMatrix< Number > & get_jump_matrix()
Definition: claw_system.C:842
std::vector< std::unique_ptr< SparseMatrix< Number > > > _boundary_condition_matrices
The "boundary condition" matrices.
Definition: claw_system.h:196
virtual void print_info()
Print out some info about the system&#39;s configuration.
Definition: claw_system.C:249
void solve_conservation_law()
Solve the conservation law.
Definition: claw_system.C:95
void set_delta_t(Real delta_t_in)
Set/get the time-step size.
Definition: claw_system.C:858
Real _LxF_constant
The constant C in the Lax-Friedrichs flux.
Definition: claw_system.h:217
void write_out_discretization_matrices()
Print discretization matrices to file in MATLAB-readable format.
Definition: claw_system.C:929
virtual ~ClawSystem()
Destructor.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static TemporalDiscretizationType string_to_enum(std::string string_type)
Convert between string and enum for TemporalDiscretizationType.
Definition: claw_system.C:67
unsigned int _write_interval
The time step interval between writing out solutions.
Definition: claw_system.h:222
unsigned int _n_time_steps
The number of time steps.
Definition: claw_system.h:212
SparseMatrix< Number > & get_mass_matrix()
Get a reference to one of the discretization matrices.
Definition: claw_system.C:816
std::unique_ptr< SparseMatrix< Number > > _mass_matrix
The mass matrix.
Definition: claw_system.h:175
void assemble_advection_matrices()
Definition: claw_system.C:335
virtual std::string system_type() const override
Definition: claw_system.C:61
void assemble_avg_coupling_matrices()
Definition: claw_system.C:415
virtual void assemble_claw_rhs(NumericVector< Number > &)=0
Assemble the right-hand side vector.
const std::string & name() const
Definition: system.h:2342
void set_LxF_constant(Real LxF_constant_in)
Set/get the Lax-Friedrichs constant.
Definition: claw_system.C:882
unsigned int get_write_interval()
Definition: claw_system.C:912
TemporalDiscretizationType _temporal_discretization_type
String that defines the type of temporal discretization that we use.
Definition: claw_system.h:202
std::vector< std::unique_ptr< SparseMatrix< Number > > > _advection_matrices
The "advection" matrices.
Definition: claw_system.h:180