The finite element method using deal.II - 2021/2022
Step4< dim > Class Template Reference

Public Member Functions

 Step4 ()
 
void run ()
 

Private Member Functions

void make_grid ()
 
void setup_system ()
 
void assemble_system ()
 
void solve ()
 
void output_results () const
 

Private Attributes

Triangulation< dim > triangulation
 
FE_Q< dim > fe
 
DoFHandler< dim > dof_handler
 
SparsityPattern sparsity_pattern
 
SparseMatrix< double > system_matrix
 
Vector< double > solution
 
Vector< double > system_rhs
 

Detailed Description

template<int dim>
class Step4< dim >

Definition at line 49 of file step-4.cc.

Constructor & Destructor Documentation

◆ Step4()

template<int dim>
Step4< dim >::Step4

Definition at line 117 of file step-4.cc.

118  : fe(1)
120 {}

Member Function Documentation

◆ run()

template<int dim>
void Step4< dim >::run

Definition at line 233 of file step-4.cc.

234 {
235  std::cout << "Solving problem in " << dim << " space dimensions."
236  << std::endl;
237  make_grid();
238  setup_system();
239  assemble_system();
240  solve();
241  output_results();
242 }

◆ make_grid()

template<int dim>
void Step4< dim >::make_grid
private

Definition at line 125 of file step-4.cc.

126 {
129  std::cout << " Number of active cells: " << triangulation.n_active_cells()
130  << std::endl
131  << " Total number of cells: " << triangulation.n_cells()
132  << std::endl;
133 }

◆ setup_system()

template<int dim>
void Step4< dim >::setup_system
private

Definition at line 138 of file step-4.cc.

139 {
141  std::cout << " Number of degrees of freedom: " << dof_handler.n_dofs()
142  << std::endl;
149 }

◆ assemble_system()

template<int dim>
void Step4< dim >::assemble_system
private

Definition at line 154 of file step-4.cc.

155 {
156  QGauss<dim> quadrature_formula(fe.degree + 1);
157  RightHandSide<dim> right_hand_side;
158  FEValues<dim> fe_values(fe,
159  quadrature_formula,
160  update_values | update_gradients |
161  update_quadrature_points | update_JxW_values);
162  const unsigned int dofs_per_cell = fe.dofs_per_cell;
163  FullMatrix<double> cell_matrix(dofs_per_cell, dofs_per_cell);
164  Vector<double> cell_rhs(dofs_per_cell);
165  std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
166  for (const auto &cell : dof_handler.active_cell_iterators())
167  {
168  fe_values.reinit(cell);
169  cell_matrix = 0;
170  cell_rhs = 0;
171  for (const unsigned int q_index : fe_values.quadrature_point_indices())
172  for (const unsigned int i : fe_values.dof_indices())
173  {
174  for (const unsigned int j : fe_values.dof_indices())
175  cell_matrix(i, j) +=
176  (fe_values.shape_grad(i, q_index) * // grad phi_i(x_q)
177  fe_values.shape_grad(j, q_index) * // grad phi_j(x_q)
178  fe_values.JxW(q_index)); // dx
179  const auto x_q = fe_values.quadrature_point(q_index);
180  cell_rhs(i) += (fe_values.shape_value(i, q_index) * // phi_i(x_q)
181  right_hand_side.value(x_q) * // f(x_q)
182  fe_values.JxW(q_index)); // dx
183  }
184  cell->get_dof_indices(local_dof_indices);
185  for (const unsigned int i : fe_values.dof_indices())
186  {
187  for (const unsigned int j : fe_values.dof_indices())
188  system_matrix.add(local_dof_indices[i],
189  local_dof_indices[j],
190  cell_matrix(i, j));
191  system_rhs(local_dof_indices[i]) += cell_rhs(i);
192  }
193  }
194  std::map<types::global_dof_index, double> boundary_values;
196  0,
198  boundary_values);
199  MatrixTools::apply_boundary_values(boundary_values,
201  solution,
202  system_rhs);
203 }

◆ solve()

template<int dim>
void Step4< dim >::solve
private

Definition at line 208 of file step-4.cc.

209 {
210  SolverControl solver_control(1000, 1e-12);
211  SolverCG<Vector<double>> solver(solver_control);
213  std::cout << " " << solver_control.last_step()
214  << " CG iterations needed to obtain convergence." << std::endl;
215 }

◆ output_results()

template<int dim>
void Step4< dim >::output_results
private

Definition at line 220 of file step-4.cc.

221 {
222  DataOut<dim> data_out;
224  data_out.add_data_vector(solution, "solution");
225  data_out.build_patches();
226  std::ofstream output(dim == 2 ? "solution-2d.vtk" : "solution-3d.vtk");
227  data_out.write_vtk(output);
228 }

Member Data Documentation

◆ triangulation

template<int dim>
Triangulation<dim> Step4< dim >::triangulation
private

Definition at line 67 of file step-4.cc.

◆ fe

template<int dim>
FE_Q<dim> Step4< dim >::fe
private

Definition at line 68 of file step-4.cc.

◆ dof_handler

template<int dim>
DoFHandler<dim> Step4< dim >::dof_handler
private

Definition at line 69 of file step-4.cc.

◆ sparsity_pattern

template<int dim>
SparsityPattern Step4< dim >::sparsity_pattern
private

Definition at line 70 of file step-4.cc.

◆ system_matrix

template<int dim>
SparseMatrix<double> Step4< dim >::system_matrix
private

Definition at line 71 of file step-4.cc.

◆ solution

template<int dim>
Vector<double> Step4< dim >::solution
private

Definition at line 72 of file step-4.cc.

◆ system_rhs

template<int dim>
Vector<double> Step4< dim >::system_rhs
private

Definition at line 73 of file step-4.cc.


The documentation for this class was generated from the following file:
SparseMatrix< double >::add
void add(const size_type i, const size_type j, const double value)
local_dof_indices
std::vector< types::global_dof_index > local_dof_indices
SolverCG
PreconditionIdentity
PreconditionIdentity()
SolverControl
RightHandSide::value
virtual double value(const Point< dim > &p, const unsigned int component=0) const override
Definition: step-4.cc:97
parallel::distributed::Triangulation::refine_global
void refine_global(const unsigned int times=1)
Step4::output_results
void output_results() const
Definition: step-4.cc:220
MatrixTools::apply_boundary_values
void apply_boundary_values(const std::map< types::global_dof_index, number > &boundary_values, SparseMatrix< number > &matrix, Vector< number > &solution, Vector< number > &right_hand_side, const bool eliminate_columns=true)
DataOut::write_vtk
void write_vtk(std::ostream &out) const
RightHandSide
Definition: step-4.cc:78
cell_matrix
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double >> &velocity, const double factor=1.)
DataOut::attach_dof_handler
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
BoundaryValues
Definition: step-4.cc:87
Step4::sparsity_pattern
SparsityPattern sparsity_pattern
Definition: step-4.cc:70
DoFHandler::distribute_dofs
void distribute_dofs(const FiniteElement< dim, spacedim > &fe)
FEValues< dim >
FE_Q::degree
const unsigned int degree
DoFTools::make_sparsity_pattern
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity_pattern, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
Step4::dof_handler
DoFHandler< dim > dof_handler
Definition: step-4.cc:69
SparseMatrix< double >::reinit
virtual void reinit(const SparsityPattern &sparsity)
FE_Q::dofs_per_cell
const unsigned int dofs_per_cell
Step4::solve
void solve()
Definition: step-4.cc:208
DataOut::build_patches
virtual void build_patches(const unsigned int n_subdivisions=0)
Step4::fe
FE_Q< dim > fe
Definition: step-4.cc:68
DynamicSparsityPattern
VectorTools::interpolate_boundary_values
void interpolate_boundary_values(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask=ComponentMask())
Step4::system_rhs
Vector< double > system_rhs
Definition: step-4.cc:73
output
output
Vector< double >::reinit
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
QGauss
Step4::assemble_system
void assemble_system()
Definition: step-4.cc:154
Step4::make_grid
void make_grid()
Definition: step-4.cc:125
SparsityPattern::copy_from
void copy_from(const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end)
GridGenerator::hyper_cube
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
DataOut::add_data_vector
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={})
Step4::triangulation
Triangulation< dim > triangulation
Definition: step-4.cc:67
Step4::solution
Vector< double > solution
Definition: step-4.cc:72
parallel::distributed::Triangulation::n_active_cells
unsigned int n_active_cells() const
Step4::setup_system
void setup_system()
Definition: step-4.cc:138
FullMatrix< double >
DoFHandler::active_cell_iterators
IteratorRange< active_cell_iterator > active_cell_iterators() const
parallel::distributed::Triangulation::n_cells
unsigned int n_cells() const
DataOut< dim >
Vector< double >
Step4::system_matrix
SparseMatrix< double > system_matrix
Definition: step-4.cc:71
DoFHandler::n_dofs
types::global_dof_index n_dofs() const