libMesh
transient_system.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2024 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 
19 
20 #ifndef LIBMESH_TRANSIENT_SYSTEM_H
21 #define LIBMESH_TRANSIENT_SYSTEM_H
22 
23 // Local Includes
24 #include "libmesh/system.h"
25 #include "libmesh/libmesh_config.h"
26 
27 namespace libMesh
28 {
29 
30 // Forward declarations
31 class LinearImplicitSystem;
32 class NonlinearImplicitSystem;
33 class ExplicitSystem;
34 #ifdef LIBMESH_HAVE_SLEPC
35 class EigenSystem;
36 #endif
37 
56 template <class Base>
57 class TransientSystem : public Base
58 {
59 public:
60 
65  const std::string & name,
66  const unsigned int number);
67 
75  TransientSystem (const TransientSystem &) = delete;
76  TransientSystem & operator= (const TransientSystem &) = delete;
77  TransientSystem (TransientSystem &&) = default;
79  virtual ~TransientSystem ();
80 
85 
89  sys_type & system () { return *this; }
90 
95  virtual void clear () override;
96 
102  virtual std::string system_type () const override;
103 
104 
105  //-----------------------------------------------------------------
106  // access to the solution data fields
107 
112  Number old_solution (const dof_id_type global_dof_number) const;
113 
118  Number older_solution (const dof_id_type global_dof_number) const;
119 
127 
135 
136 protected:
137 
143  virtual void re_update () override;
144 
145 private:
146 
150  virtual void add_old_vectors ();
151 };
152 
153 
154 
155 // -----------------------------------------------------------
156 // Useful typedefs
162 #ifdef LIBMESH_HAVE_SLEPC
164 #endif
165 
166 
167 
168 // ------------------------------------------------------------
169 // TransientSystem inline methods
170 template <class Base>
171 inline
173 {
174  std::string type = "Transient";
175  type += Base::system_type ();
176 
177  return type;
178 }
179 
180 
181 
182 } // namespace libMesh
183 
184 
185 
186 
187 #endif // LIBMESH_TRANSIENT_SYSTEM_H
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
TransientSystem< LinearImplicitSystem > TransientImplicitSystem
This is the EquationSystems class.
TransientSystem< Base > sys_type
The type of system.
TransientSystem< EigenSystem > TransientEigenSystem
virtual std::string system_type() const override
virtual void add_old_vectors()
Helper function for (re-)adding old and older solution vectors.
TransientSystem< System > TransientBaseSystem
TransientSystem< NonlinearImplicitSystem > TransientNonlinearImplicitSystem
Number older_solution(const dof_id_type global_dof_number) const
Manages storage and variables for transient systems.
The libMesh namespace provides an interface to certain functionality in the library.
NumericVector< Number > * old_local_solution
All the values I need to compute my contribution to the simulation at hand.
TransientSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
NumericVector< Number > * older_local_solution
All the values I need to compute my contribution to the simulation at hand.
TransientSystem< ExplicitSystem > TransientExplicitSystem
Number old_solution(const dof_id_type global_dof_number) const
TransientSystem< LinearImplicitSystem > TransientLinearImplicitSystem
virtual void clear() override
Clear all the data structures associated with the system.
This class is part of the rbOOmit framework.
TransientSystem & operator=(const TransientSystem &)=delete
virtual void re_update() override
Re-update the local values when the mesh has changed.
uint8_t dof_id_type
Definition: id_types.h:67