Go to the source code of this file.
◆ apply_initial()
void apply_initial |
( |
EquationSystems & |
es, |
|
|
const std::string & |
system_name |
|
) |
| |
◆ assemble_wave()
void assemble_wave |
( |
EquationSystems & |
es, |
|
|
const std::string & |
system_name |
|
) |
| |
Definition at line 323 of file transient_ex2.C.
328 libmesh_assert_equal_to (system_name,
"Wave");
346 FEType fe_type = t_system.get_dof_map().variable_type(0);
372 std::unique_ptr<FEBase> fe (FEBase::build(
dim, fe_type));
378 fe->attach_quadrature_rule (&qrule);
381 const std::vector<Real> & JxW = fe->get_JxW();
384 const std::vector<std::vector<Real>> & phi = fe->get_phi();
388 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
393 const DofMap & dof_map = t_system.get_dof_map();
403 std::vector<dof_id_type> dof_indices;
429 const unsigned int n_dof_indices = dof_indices.size();
431 Ke.
resize (n_dof_indices, n_dof_indices);
432 Ce.
resize (n_dof_indices, n_dof_indices);
433 Me.
resize (n_dof_indices, n_dof_indices);
434 zero_matrix.
resize (n_dof_indices, n_dof_indices);
435 Fe.
resize (n_dof_indices);
440 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
445 for (std::size_t i=0; i<phi.size(); i++)
446 for (std::size_t j=0; j<phi.size(); j++)
448 Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
449 Me(i,j) += JxW[qp]*phi[i][qp]*phi[j][qp]
463 for (
auto side : elem->side_index_range())
469 std::unique_ptr<FEBase> fe_face (FEBase::build(
dim, fe_type));
478 fe_face->attach_quadrature_rule (&qface);
482 const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
486 const std::vector<Real> & JxW_face = fe_face->get_JxW();
490 fe_face->reinit(elem, side);
494 const Real acc_n_value = 1.0;
497 for (
unsigned int qp=0; qp<qface.n_points(); qp++)
501 for (std::size_t i=0; i<phi_face.size(); i++)
503 Fe(i) += acc_n_value*rho
504 *phi_face[i][qp]*JxW_face[qp];
References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), dim, libMesh::DofMap::dof_indices(), libMesh::Parameters::get(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::QBase::n_points(), libMesh::EquationSystems::parameters, libMesh::Real, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), and libMesh::SECOND.
Referenced by main().
◆ fill_dirichlet_bc()
void fill_dirichlet_bc |
( |
EquationSystems & |
es, |
|
|
const std::string & |
system_name |
|
) |
| |
Definition at line 562 of file transient_ex2.C.
567 libmesh_assert_equal_to (system_name,
"Wave");
584 const bool do_for_matrix =
590 for (
unsigned int n_cnt=0; n_cnt<
n_nodes; n_cnt++)
600 const Real z_coo = 4.;
605 unsigned int dn = curr_node.
dof_number(0, 0, 0);
608 const Real penalty = 1.e10;
613 if (t_system.
time < .002)
614 p_value = sin(2*
pi*t_system.
time/.002);
619 rhs.
add(dn, p_value*penalty);
624 matrix.add(dn, dn, penalty);
References std::abs(), libMesh::NumericVector< T >::add(), libMesh::DofObject::dof_number(), libMesh::Parameters::get(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::ImplicitSystem::matrix, mesh, n_nodes, libMesh::MeshBase::n_nodes(), libMesh::MeshBase::node_ref(), libMesh::EquationSystems::parameters, libMesh::pi, libMesh::Real, libMesh::ExplicitSystem::rhs, libMesh::System::time, and libMesh::TOLERANCE.
Referenced by main().
◆ main()
int main |
( |
int |
argc, |
|
|
char ** |
argv |
|
) |
| |
Definition at line 100 of file transient_ex2.C.
107 "--enable-petsc, --enable-trilinos, or --enable-eigen");
111 libmesh_error_msg(
"Usage: " << argv[0] <<
" [meshfile]");
118 for (
int i=1; i<argc; i++)
132 std::string mesh_file = argv[1];
133 libMesh::out <<
"Mesh file is: " << mesh_file << std::endl;
136 libmesh_example_requires(3 <= LIBMESH_DIM,
"3D support");
153 const unsigned int result_node = 274;
162 const Real delta_t = .0000625;
165 unsigned int n_time_steps = 300;
192 t_system.set_newmark_parameters(delta_t);
197 equation_systems.parameters.set<
Real>(
"speed") = 1000.;
198 equation_systems.parameters.set<
Real>(
"fluid density") = 1000.;
204 equation_systems.init();
207 equation_systems.print_info();
210 std::ofstream res_out(
"pressure_node.res");
214 const unsigned int res_node_no = result_node;
216 unsigned int dof_no = res_node.
dof_number(0, 0, 0);
230 res_out <<
"# pressure at node " << res_node_no <<
"\n"
231 <<
"# time\tpressure\n"
232 << t_system.time <<
"\t" << 0 << std::endl;
235 for (
unsigned int time_step=0; time_step<n_time_steps; time_step++)
239 t_system.time += delta_t;
242 t_system.update_rhs();
257 equation_systems.parameters.set<
bool>(
"Newmark set BC for Matrix") =
true;
262 equation_systems.parameters.set<
bool>(
"Newmark set BC for Matrix") =
false;
273 if (time_step == 30 || time_step == 60 ||
274 time_step == 90 || time_step == 120)
276 std::ostringstream file_name;
278 #ifdef LIBMESH_HAVE_VTK
301 t_system.update_u_v_a();
306 std::vector<Number> global_displacement(displacement.
size());
307 displacement.
localize(global_displacement);
312 res_out << t_system.time <<
"\t"
313 << global_displacement[dof_no]
References libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), apply_initial(), assemble_wave(), libMesh::default_solver_package(), libMesh::DofObject::dof_number(), libMesh::EIGEN_SOLVERS, fill_dirichlet_bc(), libMesh::FIRST, libMesh::EquationSystems::get_system(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, libMesh::NumericVector< T >::localize(), mesh, libMesh::MeshBase::node_ref(), libMesh::out, libMesh::EquationSystems::parameters, libMesh::PETSC_SOLVERS, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::MeshBase::read(), libMesh::Real, libMesh::Parameters::set(), libMesh::NumericVector< T >::size(), and libMesh::MeshOutput< MT >::write_equation_systems().
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.
const MeshBase & get_mesh() const
void assemble_wave(EquationSystems &es, const std::string &system_name)
virtual void read(const std::string &name, void *mesh_data=nullptr, bool skip_renumber_nodes_and_elements=false, bool skip_find_neighbors=false)=0
Interfaces for reading/writing a mesh to/from a file.
NumericVector< Number > * rhs
The system matrix.
This class implements specific orders of Gauss quadrature.
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
This class contains a specific system class.
const T_sys & get_system(const std::string &name) const
static const Real TOLERANCE
SolverPackage default_solver_package()
unsigned int mesh_dimension() const
This class implements writing meshes in the GMV format.
This class implements reading and writing meshes in the VTK format.
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
virtual numeric_index_type size() const =0
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i].
This is the MeshBase class.
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.
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.
SparseMatrix< Number > * matrix
The system matrix.
A Node is like a Point, but with more information.
void apply_initial(EquationSystems &es, const std::string &system_name)
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
const dof_id_type n_nodes
This is the EquationSystems class.
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
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 ...
virtual const Node & node_ref(const dof_id_type i) const
virtual dof_id_type n_nodes() const =0
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...
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 fill_dirichlet_bc(EquationSystems &es, const std::string &system_name)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const T & get(const std::string &) const
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Parameters parameters
Data structure holding arbitrary parameters.