The finite element method using deal.II - 2021/2022
Step3 Class Reference

#include <step-3.h>

Inheritance diagram for Step3:
[legend]

Public Member Functions

 Step3 ()
 
void run ()
 

Protected Member Functions

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

Protected Attributes

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

Friends

class Step3Tester
 

Detailed Description

Definition at line 55 of file step-3.h.

Constructor & Destructor Documentation

◆ Step3()

Step3::Step3 ( )

Definition at line 24 of file step-3.cc.

25  : fe(1)
27 {}

Member Function Documentation

◆ run()

void Step3::run ( )

Definition at line 133 of file step-3.cc.

134 {
135  make_grid();
136  setup_system();
137  assemble_system();
138  solve();
139  output_results();
140 }

◆ make_grid()

void Step3::make_grid ( )
protected

Definition at line 32 of file step-3.cc.

33 {
36  std::cout << "Number of active cells: " << triangulation.n_active_cells()
37  << std::endl;
38 }

◆ setup_system()

void Step3::setup_system ( )
protected

Definition at line 43 of file step-3.cc.

44 {
46  std::cout << "Number of degrees of freedom: " << dof_handler.n_dofs()
47  << std::endl;
54 }

◆ assemble_system()

void Step3::assemble_system ( )
protected

Definition at line 59 of file step-3.cc.

60 {
61  QGauss<2> quadrature_formula(fe.degree + 1);
62  FEValues<2> fe_values(fe,
63  quadrature_formula,
64  update_values | update_gradients | update_JxW_values);
65  const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
66  FullMatrix<double> cell_matrix(dofs_per_cell, dofs_per_cell);
67  Vector<double> cell_rhs(dofs_per_cell);
68  std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
69  for (const auto &cell : dof_handler.active_cell_iterators())
70  {
71  fe_values.reinit(cell);
72  cell_matrix = 0;
73  cell_rhs = 0;
74  for (const unsigned int q_index : fe_values.quadrature_point_indices())
75  {
76  for (const unsigned int i : fe_values.dof_indices())
77  for (const unsigned int j : fe_values.dof_indices())
78  cell_matrix(i, j) +=
79  (fe_values.shape_grad(i, q_index) * // grad phi_i(x_q)
80  fe_values.shape_grad(j, q_index) * // grad phi_j(x_q)
81  fe_values.JxW(q_index)); // dx
82  for (const unsigned int i : fe_values.dof_indices())
83  cell_rhs(i) += (fe_values.shape_value(i, q_index) * // phi_i(x_q)
84  1. * // f(x_q)
85  fe_values.JxW(q_index)); // dx
86  }
87  cell->get_dof_indices(local_dof_indices);
88  for (const unsigned int i : fe_values.dof_indices())
89  for (const unsigned int j : fe_values.dof_indices())
90  system_matrix.add(local_dof_indices[i],
91  local_dof_indices[j],
92  cell_matrix(i, j));
93  for (const unsigned int i : fe_values.dof_indices())
94  system_rhs(local_dof_indices[i]) += cell_rhs(i);
95  }
96  std::map<types::global_dof_index, double> boundary_values;
98  0,
100  boundary_values);
101  MatrixTools::apply_boundary_values(boundary_values,
103  solution,
104  system_rhs);
105 }

◆ solve()

void Step3::solve ( )
protected

Definition at line 110 of file step-3.cc.

111 {
112  SolverControl solver_control(1000, 1e-12);
113  SolverCG<Vector<double>> solver(solver_control);
115 }

◆ output_results()

void Step3::output_results ( ) const
protected

Definition at line 120 of file step-3.cc.

121 {
122  DataOut<2> data_out;
124  data_out.add_data_vector(solution, "solution");
125  data_out.build_patches();
126  std::ofstream output("solution.vtk");
127  data_out.write_vtk(output);
128 }

Friends And Related Function Documentation

◆ Step3Tester

friend class Step3Tester
friend

Definition at line 82 of file step-3.h.

Member Data Documentation

◆ triangulation

Triangulation<2> Step3::triangulation
protected

Definition at line 74 of file step-3.h.

◆ fe

FE_Q<2> Step3::fe
protected

Definition at line 75 of file step-3.h.

◆ dof_handler

DoFHandler<2> Step3::dof_handler
protected

Definition at line 76 of file step-3.h.

◆ sparsity_pattern

SparsityPattern Step3::sparsity_pattern
protected

Definition at line 77 of file step-3.h.

◆ system_matrix

SparseMatrix<double> Step3::system_matrix
protected

Definition at line 78 of file step-3.h.

◆ solution

Vector<double> Step3::solution
protected

Definition at line 79 of file step-3.h.

◆ system_rhs

Vector<double> Step3::system_rhs
protected

Definition at line 80 of file step-3.h.


The documentation for this class was generated from the following files:
SparseMatrix< double >::add
void add(const size_type i, const size_type j, const double value)
Step3::setup_system
void setup_system()
Definition: step-3.cc:43
local_dof_indices
std::vector< types::global_dof_index > local_dof_indices
SolverCG
PreconditionIdentity
PreconditionIdentity()
Step3::fe
FE_Q< 2 > fe
Definition: step-3.h:75
SolverControl
Step3::system_rhs
Vector< double > system_rhs
Definition: step-3.h:80
Step3::assemble_system
void assemble_system()
Definition: step-3.cc:59
Triangulation::refine_global
void refine_global(const unsigned int times=1)
Step3::make_grid
void make_grid()
Definition: step-3.cc:32
Step3::dof_handler
DoFHandler< 2 > dof_handler
Definition: step-3.h:76
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
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 > &)
DoFHandler::distribute_dofs
void distribute_dofs(const FiniteElement< dim, spacedim > &fe)
FEValues
FE_Q::degree
const unsigned int degree
Step3::sparsity_pattern
SparsityPattern sparsity_pattern
Definition: step-3.h:77
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)
SparseMatrix< double >::reinit
virtual void reinit(const SparsityPattern &sparsity)
DataOut::build_patches
virtual void build_patches(const unsigned int n_subdivisions=0)
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())
Step3::triangulation
Triangulation< 2 > triangulation
Definition: step-3.h:74
output
output
Vector< double >::reinit
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
QGauss
Step3::output_results
void output_results() const
Definition: step-3.cc:120
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)
Step3::solve
void solve()
Definition: step-3.cc:110
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={})
Step3::system_matrix
SparseMatrix< double > system_matrix
Definition: step-3.h:78
Triangulation::n_active_cells
unsigned int n_active_cells() const
Functions::ZeroFunction
FullMatrix< double >
DoFHandler::active_cell_iterators
IteratorRange< active_cell_iterator > active_cell_iterators() const
FE_Q::n_dofs_per_cell
unsigned int n_dofs_per_cell() const
DataOut
Vector< double >
Step3::solution
Vector< double > solution
Definition: step-3.h:79
DoFHandler::n_dofs
types::global_dof_index n_dofs() const