frequency_system.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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
30 
31 // C++ includes
32 #include <string>
33 #include <vector>
34 
35 namespace libMesh
36 {
37 
64 {
65 public:
66 
72  const std::string & name_in,
73  const unsigned int number_in);
78 
85  virtual void clear () override;
86 
92  void clear_all ();
93 
98  virtual void assemble () override;
99 
103  virtual void solve () override;
104 
113  void solve (const unsigned int n_start,
114  const unsigned int n_stop);
115 
120  virtual std::string system_type () const override { return "Frequency"; }
121 
122 
123  //--------------------------------------------------------
124  // Methods specific to the FrequencySystem
125  //
126 
139  void set_frequencies_by_steps (const Number base_freq,
140  const Number freq_step=0.,
141  const unsigned int n_freq=1,
142  const bool allocate_solution_duplicates=true);
143 
153  void set_frequencies_by_range (const Number min_freq,
154  const Number max_freq,
155  const unsigned int n_freq,
156  const bool allocate_solution_duplicates=true);
157 
165  void set_frequencies (const std::vector<Real> & frequencies,
166  const bool allocate_solution_duplicates=true);
167 
168 
169  void set_frequencies (const std::vector<Number> & frequencies,
170  const bool allocate_solution_duplicates=true);
171 
175  unsigned int n_frequencies () const;
176 
183  void attach_solve_function(void fptr(EquationSystems & es,
184  const std::string & name));
185 
190  const std::string & name);
191 
195  std::pair<unsigned int, Real> get_rval (unsigned int n) const;
196 
202  std::string form_freq_param_name(const unsigned int n) const;
203 
209  std::string form_solu_vec_name(const unsigned int n) const;
210 
211 
212 protected:
213 
214 
221  virtual void init_data () override;
222 
227  void set_current_frequency(unsigned int n);
228 
235 
244 
251 
265 
270  std::vector<std::pair<unsigned int, Real>> vec_rval;
271 
272 };
273 
274 
275 
276 // ------------------------------------------------------------
277 // FrequencySystem inline methods
278 inline
279 std::pair<unsigned int, Real> FrequencySystem::get_rval (unsigned int n) const
280 {
281  libmesh_assert_less (n, vec_rval.size());
282 
283  return vec_rval[n];
284 }
285 
286 
287 } // namespace libMesh
288 
289 #endif // if defined(LIBMESH_USE_COMPLEX_NUMBERS)
290 
291 #endif // LIBMESH_FREQUENCY_SYSTEM_H
void attach_solve_function(void fptr(EquationSystems &es, const std::string &name))
virtual void solve() override
Manages multiples systems of equations.
void(* solve_system)(EquationSystems &es, const std::string &name)
virtual void assemble() override
virtual std::string system_type() const override
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
void set_frequencies_by_range(const Number min_freq, const Number max_freq, const unsigned int n_freq, const bool allocate_solution_duplicates=true)
void set_frequencies(const std::vector< Real > &frequencies, const bool allocate_solution_duplicates=true)
std::string form_freq_param_name(const unsigned int n) const
FrequencySystem(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
virtual void init_data() override
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)
std::pair< unsigned int, Real > get_rval(unsigned int n) const
virtual void clear() override
unsigned int n_frequencies() const
void set_current_frequency(unsigned int n)
const std::string & name() const
Definition: system.h:2017
std::string form_solu_vec_name(const unsigned int n) const
std::vector< std::pair< unsigned int, Real > > vec_rval