22 #ifndef _aspect_stokes_matrix_free_h 23 #define _aspect_stokes_matrix_free_h 29 #include <deal.II/matrix_free/matrix_free.h> 30 #include <deal.II/matrix_free/operators.h> 31 #include <deal.II/matrix_free/fe_evaluation.h> 33 #include <deal.II/multigrid/mg_constrained_dofs.h> 34 #include <deal.II/multigrid/multigrid.h> 35 #include <deal.II/multigrid/mg_transfer_matrix_free.h> 36 #include <deal.II/multigrid/mg_tools.h> 37 #include <deal.II/multigrid/mg_coarse.h> 38 #include <deal.II/multigrid/mg_smoother.h> 39 #include <deal.II/multigrid/mg_matrix.h> 41 #include <deal.II/lac/vector.h> 42 #include <deal.II/lac/block_vector.h> 43 #include <deal.II/lac/la_parallel_vector.h> 44 #include <deal.II/lac/la_parallel_block_vector.h> 57 namespace TangentialBoundaryFunctions
62 const types::boundary_id bid,
63 const unsigned int first_vector_component,
64 MGConstrainedDoFs &mg_constrained_dofs);
67 add_constraint(
const std::array<types::global_dof_index,dim> &dof_indices,
68 const Tensor<1, dim> &constraining_vector,
69 AffineConstraints<double> &constraints,
70 const double inhomogeneity = 0);
72 template <
int dim,
int spacedim>
74 const MGConstrainedDoFs &mg_constrained_dofs,
75 const Mapping<dim> &mapping,
76 const unsigned int level,
77 const unsigned int first_vector_component,
78 const std::set<types::boundary_id> &boundary_ids,
79 AffineConstraints<double> &constraints);
87 namespace ChangeVectorTypes
90 const ::LinearAlgebra::ReadWriteVector<double> &rwv,
91 const VectorOperation::values operation);
94 const ::LinearAlgebra::distributed::Vector<double> &in);
96 void copy(::LinearAlgebra::distributed::Vector<double> &out,
100 const ::LinearAlgebra::distributed::BlockVector<double> &in);
102 void copy(::LinearAlgebra::distributed::BlockVector<double> &out,
110 namespace MatrixFreeStokesOperators
123 template <
int dim,
typename number>
180 Table<2, SymmetricTensor<2, dim, VectorizedArray<number>>>
200 memory_consumption()
const;
211 template <
int dim,
int degree_v,
typename number>
213 :
public MatrixFreeOperators::Base<dim, ::LinearAlgebra::distributed::BlockVector<number>>
225 void clear ()
override;
237 void compute_diagonal ()
override;
245 void apply_add (::LinearAlgebra::distributed::BlockVector<number> &dst,
246 const ::LinearAlgebra::distributed::BlockVector<number> &src)
const override;
251 void local_apply (const ::MatrixFree<dim, number> &data,
252 ::LinearAlgebra::distributed::BlockVector<number> &dst,
253 const ::LinearAlgebra::distributed::BlockVector<number> &src,
254 const std::pair<unsigned int, unsigned int> &cell_range)
const;
259 void local_apply_face (const ::MatrixFree<dim, number> &data,
260 ::LinearAlgebra::distributed::BlockVector<number> &dst,
261 const ::LinearAlgebra::distributed::BlockVector<number> &src,
262 const std::pair<unsigned int, unsigned int> &face_range)
const;
267 void local_apply_boundary_face (const ::MatrixFree<dim, number> &data,
268 ::LinearAlgebra::distributed::BlockVector<number> &dst,
269 const ::LinearAlgebra::distributed::BlockVector<number> &src,
270 const std::pair<unsigned int, unsigned int> &face_range)
const;
281 template <
int dim,
int degree_p,
typename number>
283 :
public MatrixFreeOperators::Base<dim, ::LinearAlgebra::distributed::Vector<number>>
295 void clear ()
override;
307 void compute_diagonal ()
override;
315 void apply_add (::LinearAlgebra::distributed::Vector<number> &dst,
316 const ::LinearAlgebra::distributed::Vector<number> &src)
const override;
321 void local_apply (const ::MatrixFree<dim, number> &data,
322 ::LinearAlgebra::distributed::Vector<number> &dst,
323 const ::LinearAlgebra::distributed::Vector<number> &src,
324 const std::pair<unsigned int, unsigned int> &cell_range)
const;
330 void local_compute_diagonal (
const MatrixFree<dim,number> &data,
331 ::LinearAlgebra::distributed::Vector<number> &dst,
332 const unsigned int &dummy,
333 const std::pair<unsigned int,unsigned int> &cell_range)
const;
345 template <
int dim,
int degree_v,
typename number>
347 :
public MatrixFreeOperators::Base<dim, ::LinearAlgebra::distributed::Vector<number>>
359 void clear ()
override;
371 void compute_diagonal ()
override;
378 void set_diagonal (const ::LinearAlgebra::distributed::Vector<number> &diag);
386 void apply_add (::LinearAlgebra::distributed::Vector<number> &dst,
387 const ::LinearAlgebra::distributed::Vector<number> &src)
const override;
392 void local_apply (const ::MatrixFree<dim, number> &data,
393 ::LinearAlgebra::distributed::Vector<number> &dst,
394 const ::LinearAlgebra::distributed::Vector<number> &src,
395 const std::pair<unsigned int, unsigned int> &cell_range)
const;
400 void local_compute_diagonal (
const MatrixFree<dim,number> &data,
401 ::LinearAlgebra::distributed::Vector<number> &dst,
402 const unsigned int &dummy,
403 const std::pair<unsigned int,unsigned int> &cell_range)
const;
429 virtual std::pair<double,double> solve()=0;
435 virtual void setup_dofs()=0;
444 virtual void assemble()=0;
451 virtual void build_preconditioner()=0;
463 virtual const DoFHandler<dim> &
464 get_dof_handler_v ()
const = 0;
470 virtual const DoFHandler<dim> &
471 get_dof_handler_p ()
const = 0;
477 virtual const DoFHandler<dim> &
478 get_dof_handler_projection ()
const = 0;
484 virtual const AffineConstraints<double> &
485 get_constraints_v ()
const = 0;
491 virtual const AffineConstraints<double> &
492 get_constraints_p ()
const = 0;
498 virtual const MGTransferMatrixFree<dim,GMGNumberType> &
499 get_mg_transfer_A ()
const = 0;
505 virtual const MGTransferMatrixFree<dim,GMGNumberType> &
506 get_mg_transfer_S ()
const = 0;
512 virtual std::size_t get_cell_data_memory_consumption()
const = 0;
525 template<
int dim,
int velocity_degree>
546 std::pair<double,double> solve()
override;
552 void setup_dofs()
override;
561 virtual void assemble()
override;
567 void build_preconditioner()
override;
579 const DoFHandler<dim> &
580 get_dof_handler_v ()
const override;
586 const DoFHandler<dim> &
587 get_dof_handler_p ()
const override;
593 const DoFHandler<dim> &
594 get_dof_handler_projection ()
const override;
600 const AffineConstraints<double> &
601 get_constraints_v ()
const override;
607 const AffineConstraints<double> &
608 get_constraints_p ()
const override;
614 const MGTransferMatrixFree<dim,GMGNumberType> &
615 get_mg_transfer_A ()
const override;
621 const MGTransferMatrixFree<dim,GMGNumberType> &
622 get_mg_transfer_S ()
const override;
629 std::size_t get_cell_data_memory_consumption()
const override;
635 void parse_parameters (ParameterHandler &prm);
642 void evaluate_material_model();
648 void correct_stokes_rhs();
683 MGLevelObject<MatrixFreeStokesOperators::OperatorCellData<dim, GMGNumberType>>
level_cell_data;
const OperatorCellData< dim, number > * cell_data
std::set< types::boundary_id > free_surface_boundary_indicators
MGLevelObject< MatrixFreeStokesOperators::OperatorCellData< dim, GMGNumberType > > level_cell_data
DoFHandler< dim > dof_handler_p
::TrilinosWrappers::MPI::BlockVector BlockVector
MatrixFreeStokesOperators::OperatorCellData< dim, GMGNumberType > active_cell_data
MGLevelObject< GMGSchurComplementMatrixType > mg_matrices_Schur_complement
AffineConstraints< double > constraints_v
bool apply_stabilization_free_surface_faces
DoFHandler< dim > dof_handler_projection
MGLevelObject< GMGABlockMatrixType > mg_matrices_A_block
const OperatorCellData< dim, number > * cell_data
bool symmetrize_newton_system
MGConstrainedDoFs mg_constrained_dofs_A_block
void compute_no_normal_flux_constraints_box(const DoFHandler< dim > &dof, const types::boundary_id bid, const unsigned int first_vector_component, MGConstrainedDoFs &mg_constrained_dofs)
MGConstrainedDoFs mg_constrained_dofs_projection
Table< 2, Tensor< 1, dim, VectorizedArray< number > > > free_surface_stabilization_term_table
void copy(::LinearAlgebra::distributed::BlockVector< double > &out, const TrilinosWrappers::MPI::BlockVector &in)
::TrilinosWrappers::MPI::Vector Vector
MGTransferMatrixFree< dim, GMGNumberType > mg_transfer_Schur_complement
StokesMatrixType stokes_matrix
const OperatorCellData< dim, number > * cell_data
void compute_no_normal_flux_constraints_shell(const DoFHandler< dim, spacedim > &dof_handler, const MGConstrainedDoFs &mg_constrained_dofs, const Mapping< dim > &mapping, const unsigned int level, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, AffineConstraints< double > &constraints)
MGTransferMatrixFree< dim, GMGNumberType > mg_transfer_A_block
void declare_parameters(ParameterHandler &prm)
Table< 2, SymmetricTensor< 2, dim, VectorizedArray< number > > > strain_rate_table
AffineConstraints< double > constraints_p
void add_constraint(const std::array< types::global_dof_index, dim > &dof_indices, const Tensor< 1, dim > &constraining_vector, AffineConstraints< double > &constraints, const double inhomogeneity=0)
std::vector< std::shared_ptr< MatrixFree< dim, double > > > matrix_free_objects
FESystem< dim > fe_projection
Table< 2, VectorizedArray< number > > viscosity
Table< 2, SymmetricTensor< 2, dim, VectorizedArray< number > > > newton_factor_wrt_strain_rate_table
SchurComplementMatrixType Schur_complement_block_matrix
DoFHandler< dim > dof_handler_v
bool enable_newton_derivatives
Table< 2, VectorizedArray< number > > newton_factor_wrt_pressure_table
MGConstrainedDoFs mg_constrained_dofs_Schur_complement
MGLevelObject<::LinearAlgebra::distributed::Vector< GMGNumberType > > level_viscosity_vector
ABlockMatrixType A_block_matrix