20 #ifndef LIBMESH_FREQUENCY_SYSTEM_H 21 #define LIBMESH_FREQUENCY_SYSTEM_H 23 #include "libmesh/libmesh_config.h" 26 #if defined(LIBMESH_USE_COMPLEX_NUMBERS) 29 #include "libmesh/linear_implicit_system.h" 71 const std::string & name_in,
72 const unsigned int number_in);
90 virtual void clear ()
override;
108 virtual void solve ()
override;
118 void solve (
const unsigned int n_start,
119 const unsigned int n_stop);
125 virtual std::string
system_type ()
const override {
return "Frequency"; }
145 const Number freq_step=0.,
146 const unsigned int n_freq=1,
147 const bool allocate_solution_duplicates=
true);
160 const unsigned int n_freq,
161 const bool allocate_solution_duplicates=
true);
171 const bool allocate_solution_duplicates=
true);
175 const bool allocate_solution_duplicates=
true);
189 const std::string &
name));
195 const std::string &
name);
200 std::pair<unsigned int, Real>
get_rval (
unsigned int n)
const;
275 std::vector<std::pair<unsigned int, Real>>
vec_rval;
286 libmesh_assert_less (n,
vec_rval.size());
294 #endif // if defined(LIBMESH_USE_COMPLEX_NUMBERS) 296 #endif // LIBMESH_FREQUENCY_SYSTEM_H void attach_solve_function(void fptr(EquationSystems &es, const std::string &name))
Register a required user function to use in assembling/solving the system.
virtual void solve() override
Solves the system for all frequencies.
This is the EquationSystems class.
void clear_all()
The full clear method also clears the frequencies (stored as parameters of the EquationSystems object...
void(* solve_system)(EquationSystems &es, const std::string &name)
Function that computes frequency-dependent data of the system.
FrequencySystem provides a specific system class for frequency-dependent (linear) systems...
virtual void assemble() override
Assemble the linear system.
virtual std::string system_type() const override
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
bool _keep_solution_duplicates
when the solution for each frequency should be stored in an additional vector, then this bool is true...
void set_frequencies_by_range(const Number min_freq, const Number max_freq, const unsigned int n_freq, const bool allocate_solution_duplicates=true)
Set the frequency range for which the system should be solved.
void set_frequencies(const std::vector< Real > &frequencies, const bool allocate_solution_duplicates=true)
Set the frequency range by simply copying the values from frequencies.
The libMesh namespace provides an interface to certain functionality in the library.
std::string form_freq_param_name(const unsigned int n) const
FrequencySystem(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
Constructor.
Number fptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
bool _finished_assemble
true when we have finished the assemble() phase.
virtual ~FrequencySystem()
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used...
bool _finished_init
true when we have finished the init() phase.
void set_frequencies_by_steps(const Number base_freq, const Number freq_step=0., const unsigned int n_freq=1, const bool allocate_solution_duplicates=true)
Set the frequency range for which the system should be solved.
std::pair< unsigned int, Real > get_rval(unsigned int n) const
virtual void clear() override
Clear all the data structures associated with the system, but leave the frequencies untouched...
unsigned int n_frequencies() const
bool _finished_set_frequencies
true when we have frequencies to solve for.
void set_current_frequency(unsigned int n)
Sets the current frequency to the n-th entry in the vector _frequencies.
const std::string & name() const
FrequencySystem & operator=(const FrequencySystem &)=delete
std::string form_solu_vec_name(const unsigned int n) const
std::vector< std::pair< unsigned int, Real > > vec_rval
The number of iterations and the final residual when the Ax=b is solved for multiple frequencies...