22 #ifndef _aspect_stokes_matrix_free_h 23 #define _aspect_stokes_matrix_free_h 30 #include <deal.II/matrix_free/matrix_free.h> 31 #include <deal.II/matrix_free/operators.h> 32 #include <deal.II/matrix_free/fe_evaluation.h> 34 #include <deal.II/multigrid/mg_constrained_dofs.h> 35 #include <deal.II/multigrid/multigrid.h> 36 #include <deal.II/multigrid/mg_transfer_matrix_free.h> 37 #include <deal.II/multigrid/mg_transfer_global_coarsening.templates.h> 38 #include <deal.II/multigrid/mg_tools.h> 39 #include <deal.II/multigrid/mg_coarse.h> 40 #include <deal.II/multigrid/mg_smoother.h> 41 #include <deal.II/multigrid/mg_matrix.h> 43 #include <deal.II/lac/vector.h> 44 #include <deal.II/lac/block_vector.h> 45 #include <deal.II/lac/la_parallel_vector.h> 46 #include <deal.II/lac/la_parallel_block_vector.h> 62 namespace ChangeVectorTypes
65 const ::LinearAlgebra::ReadWriteVector<double> &rwv,
66 const VectorOperation::values operation);
69 const ::LinearAlgebra::distributed::Vector<double> &in);
71 void copy(::LinearAlgebra::distributed::Vector<double> &out,
75 const ::LinearAlgebra::distributed::BlockVector<double> &in);
77 void copy(::LinearAlgebra::distributed::BlockVector<double> &out,
85 namespace MatrixFreeStokesOperators
98 template <
int dim,
typename number>
155 Table<2, SymmetricTensor<2, dim, VectorizedArray<number>>>
175 memory_consumption()
const;
186 template <
int dim,
int degree_v,
typename number>
188 :
public MatrixFreeOperators::Base<dim, ::LinearAlgebra::distributed::BlockVector<number>>
200 void clear ()
override;
212 void compute_diagonal ()
override;
220 void apply_add (::LinearAlgebra::distributed::BlockVector<number> &dst,
221 const ::LinearAlgebra::distributed::BlockVector<number> &src)
const override;
226 void local_apply (const ::MatrixFree<dim, number> &data,
227 ::LinearAlgebra::distributed::BlockVector<number> &dst,
228 const ::LinearAlgebra::distributed::BlockVector<number> &src,
229 const std::pair<unsigned int, unsigned int> &cell_range)
const;
234 void local_apply_face (const ::MatrixFree<dim, number> &data,
235 ::LinearAlgebra::distributed::BlockVector<number> &dst,
236 const ::LinearAlgebra::distributed::BlockVector<number> &src,
237 const std::pair<unsigned int, unsigned int> &face_range)
const;
242 void local_apply_boundary_face (const ::MatrixFree<dim, number> &data,
243 ::LinearAlgebra::distributed::BlockVector<number> &dst,
244 const ::LinearAlgebra::distributed::BlockVector<number> &src,
245 const std::pair<unsigned int, unsigned int> &face_range)
const;
256 template <
int dim,
int degree_p,
typename number>
258 :
public MatrixFreeOperators::Base<dim, ::LinearAlgebra::distributed::Vector<number>>
270 void clear ()
override;
282 void compute_diagonal ()
override;
290 void apply_add (::LinearAlgebra::distributed::Vector<number> &dst,
291 const ::LinearAlgebra::distributed::Vector<number> &src)
const override;
296 void local_apply (const ::MatrixFree<dim, number> &data,
297 ::LinearAlgebra::distributed::Vector<number> &dst,
298 const ::LinearAlgebra::distributed::Vector<number> &src,
299 const std::pair<unsigned int, unsigned int> &cell_range)
const;
304 void inner_cell_operation(FEEvaluation<dim,
320 template <
int dim,
int degree_v,
typename number>
322 :
public MatrixFreeOperators::Base<dim, ::LinearAlgebra::distributed::Vector<number>>
334 void clear ()
override;
346 void compute_diagonal ()
override;
353 void set_diagonal (const ::LinearAlgebra::distributed::Vector<number> &diag);
360 void inner_cell_operation(FEEvaluation<dim,
364 number> &velocity)
const;
370 void cell_operation(FEEvaluation<dim,
374 number> &velocity)
const;
381 void apply_add (::LinearAlgebra::distributed::Vector<number> &dst,
382 const ::LinearAlgebra::distributed::Vector<number> &src)
const override;
387 void local_apply (const ::MatrixFree<dim, number> &data,
388 ::LinearAlgebra::distributed::Vector<number> &dst,
389 const ::LinearAlgebra::distributed::Vector<number> &src,
390 const std::pair<unsigned int, unsigned int> &cell_range)
const;
411 virtual void setup_dofs()=0;
420 virtual void assemble()=0;
427 virtual void build_preconditioner()=0;
439 virtual const DoFHandler<dim> &
440 get_dof_handler_v ()
const = 0;
446 virtual const DoFHandler<dim> &
447 get_dof_handler_p ()
const = 0;
453 virtual const DoFHandler<dim> &
454 get_dof_handler_projection ()
const = 0;
460 virtual const AffineConstraints<double> &
461 get_constraints_v ()
const = 0;
467 virtual const AffineConstraints<double> &
468 get_constraints_p ()
const = 0;
475 get_mg_transfer_A ()
const = 0;
482 get_mg_transfer_S ()
const = 0;
488 virtual std::size_t get_cell_data_memory_consumption()
const = 0;
501 template <
int dim,
int velocity_degree>
526 std::string name()
const override;
550 const bool solve_newton_system,
551 const double last_pressure_normalization_adjustment,
558 void setup_dofs()
override;
567 void assemble()
override;
573 void build_preconditioner()
override;
584 void parse_parameters (ParameterHandler &prm)
override;
590 const DoFHandler<dim> &
591 get_dof_handler_v ()
const override;
597 const DoFHandler<dim> &
598 get_dof_handler_p ()
const override;
604 const DoFHandler<dim> &
605 get_dof_handler_projection ()
const override;
611 const AffineConstraints<double> &
612 get_constraints_v ()
const override;
618 const AffineConstraints<double> &
619 get_constraints_p ()
const override;
626 get_mg_transfer_A ()
const override;
633 get_mg_transfer_S ()
const override;
640 std::size_t get_cell_data_memory_consumption()
const override;
648 void evaluate_material_model();
654 void correct_stokes_rhs();
689 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
MGConstrainedDoFs mg_constrained_dofs_projection
Table< 2, Tensor< 1, dim, VectorizedArray< number > > > free_surface_stabilization_term_table
::TrilinosWrappers::MPI::Vector Vector
StokesMatrixType stokes_matrix
const OperatorCellData< dim, number > * cell_data
::TrilinosWrappers::BlockSparseMatrix BlockSparseMatrix
void declare_parameters(ParameterHandler &prm)
Table< 2, SymmetricTensor< 2, dim, VectorizedArray< number > > > strain_rate_table
MGTransferMF< dim, GMGNumberType > mg_transfer_A_block
AffineConstraints< double > constraints_p
MGTransferMatrixFree< dim, NumberType > MGTransferMF
std::vector< std::shared_ptr< MatrixFree< dim, double > > > matrix_free_objects
MGTransferMF< dim, GMGNumberType > mg_transfer_Schur_complement
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
void copy(TrilinosWrappers::MPI::Vector &out, const ::LinearAlgebra::distributed::Vector< double > &in)
Table< 2, VectorizedArray< number > > newton_factor_wrt_pressure_table
MGConstrainedDoFs mg_constrained_dofs_Schur_complement
ABlockMatrixType A_block_matrix