Go to the source code of this file.
◆ assemble() [1/2]
void assemble |
( |
EquationSystems & |
es, |
|
|
const std::string & |
libmesh_dbg_varsystem_name |
|
) |
| |
Definition at line 131 of file amr.C.
134 libmesh_assert_equal_to (system_name,
"primary");
143 std::unique_ptr<FEBase> fe(FEBase::build(
dim, fe_type));
146 fe->attach_quadrature_rule (&qrule);
148 std::unique_ptr<FEBase> fe_face(FEBase::build(
dim, fe_type));
151 fe_face->attach_quadrature_rule(&qface);
158 const std::vector<Real> & JxW_face = fe_face->get_JxW();
159 const std::vector<Real> & JxW = fe->get_JxW();
160 const std::vector<Point> & q_point = fe->get_xyz();
161 const std::vector<std::vector<Real>> & phi = fe->get_phi();
162 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
164 std::vector<dof_id_type> dof_indices_U;
165 std::vector<dof_id_type> dof_indices_V;
173 Real vol=0., area=0.;
186 const unsigned int n_phi = cast_int<unsigned int>(phi.size());
189 Kuu.
resize (n_phi, n_phi);
191 Kvv.
resize (n_phi, n_phi);
197 for (
unsigned int gp=0; gp<qrule.n_points(); gp++)
199 for (
unsigned int i=0; i<n_phi; ++i)
207 const Real f = q_point[gp]*q_point[gp];
214 Fu(i) += JxW[gp]*f*phi[i][gp];
215 Fv(i) += JxW[gp]*f*phi[i][gp];
217 for (
unsigned int j=0; j != n_phi; ++j)
220 Kuu(i,j) += JxW[gp]*((phi[i][gp])*(phi[j][gp]));
222 Kvv(i,j) += JxW[gp]*((phi[i][gp])*(phi[j][gp]) +
223 1.*((dphi[i][gp])*(dphi[j][gp])));
233 for (
auto side : elem->side_index_range())
234 if (elem->neighbor_ptr(side) ==
nullptr)
236 fe_face->reinit (elem, side);
238 for (
const auto & val : JxW_face)
References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), libMesh::DofMap::constrain_element_matrix_and_vector(), dim, libMesh::DofMap::dof_indices(), libMesh::FIFTH, libMesh::FIRST, libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::ImplicitSystem::matrix, mesh, libMesh::MeshBase::mesh_dimension(), libMesh::QBase::n_points(), libMesh::out, libMesh::Real, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), and libMesh::ExplicitSystem::rhs.
◆ assemble() [2/2]
Definition at line 233 of file miscellaneous_ex4.C.
239 #ifdef LIBMESH_ENABLE_AMR
242 libmesh_assert_equal_to (system_name,
"System");
262 std::unique_ptr<FEBase> fe (FEBase::build(
dim, fe_type));
263 std::unique_ptr<FEBase> fe_face (FEBase::build(
dim, fe_type));
267 QGauss qrule (
dim, fe_type.default_quadrature_order());
268 QGauss qface (
dim-1, fe_type.default_quadrature_order());
271 fe->attach_quadrature_rule (&qrule);
272 fe_face->attach_quadrature_rule (&qface);
277 const std::vector<Real> & JxW = fe->get_JxW();
278 const std::vector<Real> & JxW_face = fe_face->get_JxW();
281 const std::vector<std::vector<Real>> & phi = fe->get_phi();
282 const std::vector<std::vector<Real>> & psi = fe_face->get_phi();
286 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
312 std::vector<dof_id_type> dof_indices;
339 const unsigned int n_dofs =
340 cast_int<unsigned int>(dof_indices.size());
342 Ke.
resize (n_dofs, n_dofs);
354 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
357 for (
unsigned int i=0; i<n_dofs; i++)
360 Fe(i) += JxW[qp]*phi[i][qp];
362 for (
unsigned int j=0; j<n_dofs; j++)
367 (dphi[i][qp]*dphi[j][qp])
372 Ve(i) += JxW[qp]*phi[i][qp];
373 We(i) += JxW[qp]*phi[i][qp];
388 const Real penalty = 1.e10;
393 for (
auto s : elem->side_index_range())
394 if (elem->neighbor_ptr(s) ==
nullptr)
396 fe_face->reinit(elem, s);
398 for (
unsigned int qp=0; qp<qface.n_points(); qp++)
401 for (
unsigned int i=0; i<n_dofs; i++)
402 for (
unsigned int j=0; j<n_dofs; j++)
403 Ke(i,j) += penalty*JxW_face[qp]*psi[i][qp]*psi[j][qp];
423 std::vector<dof_id_type> dof_indices_backup(dof_indices);
425 dof_indices = dof_indices_backup;
449 #endif // #ifdef LIBMESH_ENABLE_AMR
References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), libMesh::SparseMatrix< T >::close(), libMesh::NumericVector< T >::close(), libMesh::DofMap::constrain_element_dyad_matrix(), libMesh::DofMap::constrain_element_matrix_and_vector(), dim, libMesh::DofMap::dof_indices(), libMesh::System::get_dof_map(), libMesh::ImplicitSystem::get_matrix(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::System::get_vector(), libMesh::libmesh_ignore(), libMesh::ImplicitSystem::matrix, mesh, libMesh::MeshBase::mesh_dimension(), libMesh::Real, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::ExplicitSystem::rhs, and libMesh::System::variable_type().
Referenced by main().
◆ main()
int main |
( |
int |
argc, |
|
|
char ** |
argv |
|
) |
| |
Definition at line 46 of file amr.C.
51 libmesh_error_msg(
"Usage: ./prog -d DIM filename");
54 const unsigned char dim = cast_int<unsigned char>(atoi(argv[2]));
56 std::string meshname (argv[3]);
70 mesh_refinement.refine_and_coarsen_elements ();
71 mesh_refinement.uniformly_refine (2);
105 mesh_refinement.uniformly_refine (1);
References libMesh::DofMap::_dof_coupling, libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), assemble(), libMesh::System::attach_assemble_function(), dim, libMesh::MeshBase::elem_ref(), libMesh::FIRST, libMesh::System::get_dof_map(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), mesh, libMesh::DofMap::print_dof_constraints(), libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::MeshBase::read(), libMesh::Elem::REFINE, libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::EquationSystems::reinit(), libMesh::CouplingMatrix::resize(), libMesh::Elem::set_refinement_flag(), libMesh::LinearImplicitSystem::solve(), libMesh::MeshRefinement::uniformly_refine(), libMesh::GMVIO::write(), and libMesh::MeshOutput< MT >::write_equation_systems().
virtual void close()=0
Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across p...
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
const MeshBase & get_mesh() const
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.
void print_dof_constraints(std::ostream &os=libMesh::out, bool print_nonlocal=false) const
Prints (from processor 0) all DoF and Node constraints.
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
void resize(const unsigned int n)
Resizes the matrix and initializes all entries to be 0.
void assemble(EquationSystems &es, const std::string &system_name)
virtual const Elem & elem_ref(const dof_id_type i) const
virtual void close()=0
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
const T_sys & get_system(const std::string &name) const
virtual void solve() override
Assembles & solves the linear system A*x=b.
unsigned int mesh_dimension() const
This class implements writing meshes in the GMV format.
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.
Implements (adaptive) mesh refinement algorithms for a MeshBase.
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.
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.
void constrain_element_dyad_matrix(DenseVector< Number > &v, DenseVector< Number > &w, std::vector< dof_id_type > &row_dofs, bool asymmetric_constraint_rows=true) const
Constrains a dyadic element matrix B = v w'.
void libmesh_ignore(const Args &...)
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.
virtual void write(const std::string &) override
This method implements writing a mesh to a specified file.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
This is the EquationSystems class.
void constrain_element_matrix_and_vector(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true) const
Constrains the element matrix and vector.
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 ...
const FEType & variable_type(const unsigned int i) const
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.
const DofMap & get_dof_map() const
void set_refinement_flag(const RefinementState rflag)
Sets the value of the refinement flag for the element.
CouplingMatrix * _dof_coupling
Degree of freedom coupling.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
const SparseMatrix< Number > & get_matrix(const std::string &mat_name) const
const NumericVector< Number > & get_vector(const std::string &vec_name) const