Go to the documentation of this file.
   45 #include "libmesh/libmesh.h" 
   46 #include "libmesh/libmesh_logging.h" 
   47 #include "libmesh/mesh.h" 
   48 #include "libmesh/mesh_generation.h" 
   49 #include "libmesh/exodusII_io.h" 
   50 #include "libmesh/unv_io.h" 
   51 #include "libmesh/equation_systems.h" 
   52 #include "libmesh/elem.h" 
   53 #include "libmesh/enum_xdr_mode.h" 
   57 #include "libmesh/frequency_system.h" 
   60 #include "libmesh/fe.h" 
   63 #include "libmesh/quadrature_gauss.h" 
   66 #include "libmesh/dense_matrix.h" 
   67 #include "libmesh/dense_vector.h" 
   70 #include "libmesh/sparse_matrix.h" 
   71 #include "libmesh/numeric_vector.h" 
   74 #include "libmesh/dof_map.h" 
   81 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 
   87                         const std::string & system_name);
 
   93                          const std::string & system_name);
 
  103 int main (
int argc, 
char ** argv)
 
  110 #ifndef LIBMESH_USE_COMPLEX_NUMBERS 
  111   libmesh_example_requires(
false, 
"--enable-complex");
 
  115   libmesh_example_requires(2 <= LIBMESH_DIM, 
"2D support");
 
  119     libmesh_error_msg(
"Usage: " << argv[0] << 
" -f [real_frequency imag_frequency]");
 
  121   if (
init.comm().size() > 1)
 
  123       if (
init.comm().rank() == 0)
 
  125           libMesh::err << 
"TODO: This example should be able to run in parallel." 
  136       for (
int i=1; i<argc; i++)
 
  145   const Number frequency_in (atof(argv[2]), atof(argv[3]));
 
  150   const unsigned int n_frequencies = 3;
 
  159   unvio.read(
"lshape.unv");
 
  172   unvio.read_dataset(
"lshape_data.unv");
 
  222   equation_systems.
init ();
 
  234         const std::vector<Number> * nodal_data = unvio.get_data(node);
 
  239             unsigned int dn = node->dof_number(0,
 
  242             freq_indep_rhs.
set(dn, (*nodal_data)[0]);
 
  250   for (
unsigned int n=0; n < n_frequencies; n++)
 
  258       f_system.
solve (n, n);
 
  262 #ifdef LIBMESH_HAVE_EXODUS_API 
  263       std::ostringstream file_name;
 
  280   equation_systems.
write (
"eqn_sys.dat", 
WRITE);
 
  289 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 
  294                         const std::string & system_name)
 
  299   libmesh_assert_equal_to (system_name, 
"Helmholtz");
 
  318   const FEType & fe_type = dof_map.variable_type(0);
 
  341   fe->attach_quadrature_rule (&qrule);
 
  344   const std::vector<Real> & JxW = fe->get_JxW();
 
  347   const std::vector<std::vector<Real>> & phi = fe->get_phi();
 
  351   const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
 
  368   std::vector<dof_id_type> dof_indices;
 
  375       START_LOG(
"elem init", 
"assemble_helmholtz");
 
  381       dof_map.dof_indices (elem, dof_indices);
 
  393         const unsigned int n_dof_indices = dof_indices.size();
 
  395         Ke.
resize          (n_dof_indices, n_dof_indices);
 
  396         Ce.
resize          (n_dof_indices, n_dof_indices);
 
  397         Me.
resize          (n_dof_indices, n_dof_indices);
 
  398         zero_matrix.
resize (n_dof_indices, n_dof_indices);
 
  399         Fe.
resize          (n_dof_indices);
 
  403       STOP_LOG(
"elem init", 
"assemble_helmholtz");
 
  407       START_LOG(
"stiffness & mass", 
"assemble_helmholtz");
 
  409       for (
unsigned int qp=0; qp<qrule.
n_points(); qp++)
 
  414           for (std::size_t i=0; i<phi.size(); i++)
 
  415             for (std::size_t j=0; j<phi.size(); j++)
 
  417                 Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
 
  418                 Me(i,j) += JxW[qp]*(phi[i][qp]*phi[j][qp]);
 
  422       STOP_LOG(
"stiffness & mass", 
"assemble_helmholtz");
 
  431       for (
auto side : elem->side_index_range())
 
  432         if (elem->neighbor_ptr(side) == 
nullptr)
 
  434             LOG_SCOPE(
"damping", 
"assemble_helmholtz");
 
  447             fe_face->attach_quadrature_rule (&qface);
 
  451             const std::vector<std::vector<Real>> & phi_face =
 
  456             const std::vector<Real> & JxW_face = fe_face->get_JxW();
 
  460             fe_face->reinit(elem, side);
 
  464             const Real an_value = 1.;
 
  467             for (
unsigned int qp=0; qp<qface.
n_points(); qp++)
 
  471                 for (std::size_t i=0; i<phi_face.size(); i++)
 
  472                   for (std::size_t j=0; j<phi_face.size(); j++)
 
  473                     Ce(i,j) += rho * an_value * JxW_face[qp] * phi_face[i][qp] * phi_face[j][qp];
 
  488       stiffness.add_matrix (Ke, dof_indices);
 
  503                          const std::string & system_name)
 
  505   START_LOG(
"init phase", 
"add_M_C_K_helmholtz");
 
  508   libmesh_assert_equal_to (system_name, 
"Helmholtz");
 
  521   const Number k     = omega / speed;
 
  540   const Number scale_stiffness (1., 0.);
 
  541   const Number scale_damping=unit_i*omega; 
 
  542   const Number scale_mass=-k*k;
 
  543   const Number scale_rhs=-unit_i*rho*omega;
 
  559   STOP_LOG(
"init phase", 
"add_M_C_K_helmholtz");
 
  561   START_LOG(
"global matrix & vector additions", 
"add_M_C_K_helmholtz");
 
  572   matrix.
add (scale_stiffness, stiffness);
 
  573   matrix.
add (scale_mass,      mass);
 
  574   matrix.
add (scale_damping,   damping);
 
  575   rhs.
add    (scale_rhs,       freq_indep_rhs);
 
  577   STOP_LOG(
"global matrix & vector additions", 
"add_M_C_K_helmholtz");
 
  582 #endif // LIBMESH_USE_COMPLEX_NUMBERS 
  
virtual void close()=0
Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across p...
 
virtual void zero()=0
Set all entries to zero.
 
virtual void add(const numeric_index_type i, const T value)=0
Adds value to each entry of the vector.
 
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
 
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
 
const MeshBase & get_mesh() const
 
NumericVector< Number > * rhs
The system matrix.
 
This class implements specific orders of Gauss quadrature.
 
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 SimpleRange< element_iterator > active_local_element_ptr_range()=0
 
unsigned int n_points() const
 
The libMesh namespace provides an interface to certain functionality in the library.
 
virtual void close()=0
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
 
SparseMatrix< Number > & add_matrix(const std::string &mat_name)
Adds the additional matrix mat_name to this system.
 
const T_sys & get_system(const std::string &name) const
 
int main(int argc, char **argv)
 
unsigned int mesh_dimension() const
 
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
 
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
 
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
 
void attach_assemble_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in assembling the system matrix and RHS.
 
This is the MeshBase class.
 
The UNVIO class implements the Ideas UNV universal file format.
 
virtual SimpleRange< node_iterator > node_ptr_range()=0
 
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
 
unsigned int add_variable(const std::string &var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
 
virtual void init()
Initialize all the systems.
 
SparseMatrix< Number > * matrix
The system matrix.
 
FrequencySystem provides a specific system class for frequency-dependent (linear) systems.
 
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
 
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value)=0
Add value to the element (i,j).
 
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
 
This is the EquationSystems class.
 
T & set(const std::string &)
 
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
 
void resize(const unsigned int n)
Resize the vector.
 
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
Add the full matrix dm to the SparseMatrix.
 
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
 
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.
 
virtual void set(const numeric_index_type i, const T value)=0
Sets v(i) = value.
 
This class handles the numbering of degrees of freedom on a mesh.
 
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
 
void assemble_helmholtz(EquationSystems &es, const std::string &system_name)
 
void write(const std::string &name, const XdrMODE, const unsigned int write_flags=(WRITE_DATA), bool partition_agnostic=true) const
Write the systems to disk using the XDR data format.
 
virtual void solve() override
Solves the system for all frequencies.
 
void add_M_C_K_helmholtz(EquationSystems &es, const std::string &system_name)
 
NumericVector< Number > & add_vector(const std::string &vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Adds the additional vector vec_name to this system.
 
const DofMap & get_dof_map() const
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
void prepare_for_use(const bool skip_renumber_nodes_and_elements=false, const bool skip_find_neighbors=false)
Prepare a newly ecreated (or read) mesh for use.
 
const T & get(const std::string &) const
 
const SparseMatrix< Number > & get_matrix(const std::string &mat_name) const
 
const NumericVector< Number > & get_vector(const std::string &vec_name) const
 
virtual void zero()=0
Set all entries to 0.
 
Parameters parameters
Data structure holding arbitrary parameters.