libMesh
frequency_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 
19 
20 #ifndef LIBMESH_FREQUENCY_SYSTEM_H
21 #define LIBMESH_FREQUENCY_SYSTEM_H
22 
23 #include "libmesh/libmesh_config.h"
24 
25 // Frequency domain solutions only possible with complex arithmetic
26 #if defined(LIBMESH_USE_COMPLEX_NUMBERS)
27 
28 // Local Includes
29 #include "libmesh/linear_implicit_system.h"
30 
31 // C++ includes
32 #include <string>
33 #include <vector>
34 
35 namespace libMesh
36 {
37 
64 {
65 public:
66 
71  const std::string & name_in,
72  const unsigned int number_in);
78  FrequencySystem (const FrequencySystem &) = delete;
79  FrequencySystem & operator= (const FrequencySystem &) = delete;
80  FrequencySystem (FrequencySystem &&) = default;
82  virtual ~FrequencySystem ();
83 
90  virtual void clear () override;
91 
97  void clear_all ();
98 
103  virtual void assemble () override;
104 
108  virtual void solve () override;
109 
118  void solve (const unsigned int n_start,
119  const unsigned int n_stop);
120 
125  virtual std::string system_type () const override { return "Frequency"; }
126 
127 
128  //--------------------------------------------------------
129  // Methods specific to the FrequencySystem
130  //
131 
144  void set_frequencies_by_steps (const Number base_freq,
145  const Number freq_step=0.,
146  const unsigned int n_freq=1,
147  const bool allocate_solution_duplicates=true);
148 
158  void set_frequencies_by_range (const Number min_freq,
159  const Number max_freq,
160  const unsigned int n_freq,
161  const bool allocate_solution_duplicates=true);
162 
170  void set_frequencies (const std::vector<Real> & frequencies,
171  const bool allocate_solution_duplicates=true);
172 
173 
174  void set_frequencies (const std::vector<Number> & frequencies,
175  const bool allocate_solution_duplicates=true);
176 
180  unsigned int n_frequencies () const;
181 
189  const std::string & name));
190 
195  const std::string & name);
196 
200  std::pair<unsigned int, Real> get_rval (unsigned int n) const;
201 
207  std::string form_freq_param_name(const unsigned int n) const;
208 
214  std::string form_solu_vec_name(const unsigned int n) const;
215 
216 
217 protected:
218 
219 
226  virtual void init_data () override;
227 
232  void set_current_frequency(unsigned int n);
233 
240 
249 
256 
270 
275  std::vector<std::pair<unsigned int, Real>> vec_rval;
276 
277 };
278 
279 
280 
281 // ------------------------------------------------------------
282 // FrequencySystem inline methods
283 inline
284 std::pair<unsigned int, Real> FrequencySystem::get_rval (unsigned int n) const
285 {
286  libmesh_assert_less (n, vec_rval.size());
287 
288  return vec_rval[n];
289 }
290 
291 
292 } // namespace libMesh
293 
294 #endif // if defined(LIBMESH_USE_COMPLEX_NUMBERS)
295 
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)
Definition: projection.C:80
bool _finished_assemble
true when we have finished the assemble() phase.
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
Definition: system.h:2342
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...