MAST
MAST::PKFlutterSolver Class Reference

#include <pk_flutter_solver.h>

Inheritance diagram for MAST::PKFlutterSolver:
Collaboration diagram for MAST::PKFlutterSolver:

Public Member Functions

 PKFlutterSolver ()
 
virtual ~PKFlutterSolver ()
 
virtual void clear ()
 clears the solution and other data from this solver More...
 
virtual void clear_solutions ()
 clears the solutions stored from a previous analysis. More...
 
void initialize (MAST::Parameter &V_param, MAST::Parameter &kr_param, MAST::Parameter &bref_param, Real rho, Real V_lower, Real V_upper, unsigned int n_V_divs, Real kr_lower, Real kr_upper, unsigned int n_kr_divs, std::vector< libMesh::NumericVector< Real > * > &basis)
 initializes the data structres for a flutter solution. More...
 
virtual unsigned int n_roots_found () const
 finds the number of critical points already identified in the procedure. More...
 
const MAST::FlutterRootBaseget_root (const unsigned int n) const
 returns the n th root in terms of ascending velocity that is found by the solver More...
 
virtual std::pair< bool, MAST::FlutterRootBase * > find_next_root (const Real g_tol, const unsigned int n_bisection_iters)
 Looks through the list of flutter cross-over points and iteratively zooms in to find the cross-over point. More...
 
virtual std::pair< bool, MAST::FlutterRootBase * > find_critical_root (const Real g_tol, const unsigned int n_bisection_iters)
 This method checks if the flutter root corresponding to the lowest velocity crossover has been calculated. More...
 
virtual void print_sorted_roots ()
 Prints the sorted roots to the output. More...
 
virtual void print_crossover_points ()
 Prints the crossover points output. More...
 
virtual void scan_for_roots ()
 Scans for flutter roots in the analyzed points, and identified the divergence (if k_red = 0. More...
 
virtual void calculate_sensitivity (MAST::FlutterRootBase &root, const MAST::FunctionBase &f)
 Calculate the sensitivity of the flutter root with respect to the parameter f. More...
 
- Public Member Functions inherited from MAST::FlutterSolverBase
 FlutterSolverBase ()
 defalut constructor More...
 
virtual ~FlutterSolverBase ()
 
void attach_assembly (MAST::StructuralFluidInteractionAssembly &assembly)
 attaches the assembly object to this solver. More...
 
void attach_steady_solver (MAST::FlutterSolverBase::SteadySolver &solver)
 attaches the steady solution object More...
 
virtual void clear_assembly_object ()
 clears the assembly object More...
 
void initialize (std::vector< libMesh::NumericVector< Real > * > &basis)
 initializes the data structres for a flutter solution. More...
 
void set_output_file (const std::string &nm)
 

Protected Member Functions

void _insert_new_solution (const Real k_red_ref, MAST::FlutterSolutionBase *sol)
 
virtual std::pair< bool, MAST::FlutterSolutionBase * > _bisection_search (const std::pair< MAST::FlutterSolutionBase *, MAST::FlutterSolutionBase * > &ref_sol_range, const unsigned int root_num, const Real g_tol, const unsigned int max_iters)
 
virtual std::unique_ptr< MAST::FlutterSolutionBase_analyze (const Real k_red, const Real v_ref, const MAST::FlutterSolutionBase *prev_sol=nullptr)
 Newton method to look for cross-over point method search. More...
 
void _initialize_matrices (const Real k_red, const Real v_ref, ComplexMatrixX &L, ComplexMatrixX &R, RealMatrixX &stiff)
 initializes the matrices for the specified k_red. More...
 
void _initialize_matrix_sensitivity_for_param (const MAST::FunctionBase &f, const Real k_red, const Real U_inf, ComplexMatrixX &L, ComplexMatrixX &R)
 Assembles the reduced order system structural and aerodynmaic matrices for specified flight velocity U_inf. More...
 
virtual void _identify_crossover_points ()
 identifies all cross-over and divergence points from analyzed roots More...
 

Protected Attributes

MAST::Parameter_velocity_param
 Parameter that define the velocity. More...
 
MAST::Parameter_kred_param
 Parameter that define the reduced frequency. More...
 
Real _rho
 flight density More...
 
MAST::Parameter_bref_param
 reference chord More...
 
std::pair< Real, Real_V_range
 range of reference values within which to find flutter roots More...
 
unsigned int _n_V_divs
 number of division in the reference value range for initial scanning More...
 
std::pair< Real, Real_kr_range
 range of reference values within which to find flutter roots More...
 
unsigned int _n_k_red_divs
 number of division in the reference value range for initial scanning More...
 
std::map< Real, MAST::FlutterSolutionBase * > _flutter_solutions
 map of velocity sorted flutter solutions More...
 
std::multimap< Real, MAST::FlutterRootCrossoverBase * > _flutter_crossovers
 the map of flutter crossover points versus average velocity of the two bounding roots More...
 
- Protected Attributes inherited from MAST::FlutterSolverBase
MAST::StructuralFluidInteractionAssembly_assembly
 structural assembly that provides the assembly of the system matrices. More...
 
std::vector< libMesh::NumericVector< Real > * > * _basis_vectors
 basis vector used to define the reduced order model More...
 
std::ofstream * _output
 file to which the result will be written More...
 
MAST::FlutterSolverBase::SteadySolver_steady_solver
 object provides the steady state solution. More...
 

Detailed Description

Definition at line 36 of file pk_flutter_solver.h.

Constructor & Destructor Documentation

MAST::PKFlutterSolver::PKFlutterSolver ( )

Definition at line 31 of file pk_flutter_solver.cpp.

MAST::PKFlutterSolver::~PKFlutterSolver ( )
virtual

Definition at line 44 of file pk_flutter_solver.cpp.

Here is the call graph for this function:

Member Function Documentation

std::unique_ptr< MAST::FlutterSolutionBase > MAST::PKFlutterSolver::_analyze ( const Real  k_red,
const Real  v_ref,
const MAST::FlutterSolutionBase prev_sol = nullptr 
)
protectedvirtual

Newton method to look for cross-over point method search.

performs an eigensolution at the specified reduced frequency, and sort the roots based on the provided solution pointer. If the pointer is nullptr, then no sorting is performed solve the eigenproblem

\[ L x = \lambda R x \]

Definition at line 839 of file pk_flutter_solver.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

std::pair< bool, MAST::FlutterSolutionBase * > MAST::PKFlutterSolver::_bisection_search ( const std::pair< MAST::FlutterSolutionBase *, MAST::FlutterSolutionBase * > &  ref_sol_range,
const unsigned int  root_num,
const Real  g_tol,
const unsigned int  max_iters 
)
protectedvirtual

Definition at line 311 of file pk_flutter_solver.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::PKFlutterSolver::_identify_crossover_points ( )
protectedvirtual

identifies all cross-over and divergence points from analyzed roots

Definition at line 711 of file pk_flutter_solver.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::PKFlutterSolver::_initialize_matrices ( const Real  k_red,
const Real  v_ref,
ComplexMatrixX L,
ComplexMatrixX R,
RealMatrixX stiff 
)
protected

initializes the matrices for the specified k_red.

UG does not account for structural damping.

Definition at line 969 of file pk_flutter_solver.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::PKFlutterSolver::_initialize_matrix_sensitivity_for_param ( const MAST::FunctionBase f,
const Real  k_red,
const Real  U_inf,
ComplexMatrixX L,
ComplexMatrixX R 
)
protected

Assembles the reduced order system structural and aerodynmaic matrices for specified flight velocity U_inf.

Definition at line 1029 of file pk_flutter_solver.cpp.

Here is the caller graph for this function:

void MAST::PKFlutterSolver::_insert_new_solution ( const Real  k_red_ref,
MAST::FlutterSolutionBase sol 
)
protected

Definition at line 662 of file pk_flutter_solver.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::PKFlutterSolver::calculate_sensitivity ( MAST::FlutterRootBase root,
const MAST::FunctionBase f 
)
virtual

Calculate the sensitivity of the flutter root with respect to the parameter f.

Definition at line 877 of file pk_flutter_solver.cpp.

void MAST::PKFlutterSolver::clear ( )
virtual

clears the solution and other data from this solver

Reimplemented from MAST::FlutterSolverBase.

Definition at line 52 of file pk_flutter_solver.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::PKFlutterSolver::clear_solutions ( )
virtual

clears the solutions stored from a previous analysis.

Definition at line 103 of file pk_flutter_solver.cpp.

Here is the caller graph for this function:

std::pair< bool, MAST::FlutterRootBase * > MAST::PKFlutterSolver::find_critical_root ( const Real  g_tol,
const unsigned int  n_bisection_iters 
)
virtual

This method checks if the flutter root corresponding to the lowest velocity crossover has been calculated.

If not, then it attempts to find that root using an iterative approach

Definition at line 457 of file pk_flutter_solver.cpp.

Here is the call graph for this function:

std::pair< bool, MAST::FlutterRootBase * > MAST::PKFlutterSolver::find_next_root ( const Real  g_tol,
const unsigned int  n_bisection_iters 
)
virtual

Looks through the list of flutter cross-over points and iteratively zooms in to find the cross-over point.

This should be called only after scan_for_roots() has been called. Potential cross-over points are sorted with increasing velocity, and this method will attempt to identify the next critical root in the order.

Definition at line 412 of file pk_flutter_solver.cpp.

Here is the call graph for this function:

const MAST::FlutterRootBase & MAST::PKFlutterSolver::get_root ( const unsigned int  n) const

returns the n th root in terms of ascending velocity that is found by the solver

Definition at line 388 of file pk_flutter_solver.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::PKFlutterSolver::initialize ( MAST::Parameter V_param,
MAST::Parameter kr_param,
MAST::Parameter bref_param,
Real  rho,
Real  V_lower,
Real  V_upper,
unsigned int  n_V_divs,
Real  kr_lower,
Real  kr_upper,
unsigned int  n_kr_divs,
std::vector< libMesh::NumericVector< Real > * > &  basis 
)

initializes the data structres for a flutter solution.

Definition at line 71 of file pk_flutter_solver.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

unsigned int MAST::PKFlutterSolver::n_roots_found ( ) const
virtual

finds the number of critical points already identified in the procedure.

Definition at line 125 of file pk_flutter_solver.cpp.

Here is the caller graph for this function:

void MAST::PKFlutterSolver::print_crossover_points ( )
virtual

Prints the crossover points output.

If no pointer to output is given then the output defined by set_output_file() is used.

Implements MAST::FlutterSolverBase.

Definition at line 285 of file pk_flutter_solver.cpp.

void MAST::PKFlutterSolver::print_sorted_roots ( )
virtual

Prints the sorted roots to the output.

Implements MAST::FlutterSolverBase.

Definition at line 213 of file pk_flutter_solver.cpp.

Here is the call graph for this function:

void MAST::PKFlutterSolver::scan_for_roots ( )
virtual

Scans for flutter roots in the analyzed points, and identified the divergence (if k_red = 0.

is specified) and flutter crossover points. The roots are organized in terms of increasing velocity.

Definition at line 142 of file pk_flutter_solver.cpp.

Here is the call graph for this function:

Member Data Documentation

MAST::Parameter* MAST::PKFlutterSolver::_bref_param
protected

reference chord

Definition at line 223 of file pk_flutter_solver.h.

std::multimap<Real, MAST::FlutterRootCrossoverBase*> MAST::PKFlutterSolver::_flutter_crossovers
protected

the map of flutter crossover points versus average velocity of the two bounding roots

Definition at line 256 of file pk_flutter_solver.h.

std::map<Real, MAST::FlutterSolutionBase*> MAST::PKFlutterSolver::_flutter_solutions
protected

map of velocity sorted flutter solutions

Definition at line 250 of file pk_flutter_solver.h.

std::pair<Real, Real> MAST::PKFlutterSolver::_kr_range
protected

range of reference values within which to find flutter roots

Definition at line 239 of file pk_flutter_solver.h.

MAST::Parameter* MAST::PKFlutterSolver::_kred_param
protected

Parameter that define the reduced frequency.

Definition at line 213 of file pk_flutter_solver.h.

unsigned int MAST::PKFlutterSolver::_n_k_red_divs
protected

number of division in the reference value range for initial scanning

Definition at line 245 of file pk_flutter_solver.h.

unsigned int MAST::PKFlutterSolver::_n_V_divs
protected

number of division in the reference value range for initial scanning

Definition at line 234 of file pk_flutter_solver.h.

Real MAST::PKFlutterSolver::_rho
protected

flight density

Definition at line 218 of file pk_flutter_solver.h.

std::pair<Real, Real> MAST::PKFlutterSolver::_V_range
protected

range of reference values within which to find flutter roots

Definition at line 228 of file pk_flutter_solver.h.

MAST::Parameter* MAST::PKFlutterSolver::_velocity_param
protected

Parameter that define the velocity.

Definition at line 208 of file pk_flutter_solver.h.


The documentation for this class was generated from the following files: