ExaDG
Loading...
Searching...
No Matches
operator_base.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef EXADG_OPERATORS_OPERATOR_BASE_H_
23#define EXADG_OPERATORS_OPERATOR_BASE_H_
24
25// deal.II
26#include <deal.II/base/subscriptor.h>
27#include <deal.II/dofs/dof_handler.h>
28#include <deal.II/lac/affine_constraints.h>
29#include <deal.II/lac/la_parallel_vector.h>
30#include <deal.II/lac/lapack_full_matrix.h>
31#ifdef DEAL_II_WITH_TRILINOS
32# include <deal.II/lac/trilinos_sparse_matrix.h>
33#endif
34#ifdef DEAL_II_WITH_PETSC
35# include <deal.II/lac/petsc_sparse_matrix.h>
36#endif
37#include <deal.II/matrix_free/matrix_free.h>
38
39// ExaDG
40#include <exadg/matrix_free/categorization.h>
41#include <exadg/matrix_free/integrators.h>
42
43#include <exadg/solvers_and_preconditioners/preconditioners/elementwise_preconditioners.h>
44#include <exadg/solvers_and_preconditioners/preconditioners/enum_types.h>
45#include <exadg/solvers_and_preconditioners/solvers/enum_types.h>
46#include <exadg/solvers_and_preconditioners/solvers/wrapper_elementwise_solvers.h>
47#include <exadg/solvers_and_preconditioners/utilities/invert_diagonal.h>
48
49#include <exadg/utilities/lazy_ptr.h>
50
51#include <exadg/operators/elementwise_operator.h>
52#include <exadg/operators/enum_types.h>
53#include <exadg/operators/integrator_flags.h>
54#include <exadg/operators/mapping_flags.h>
55#include <exadg/operators/operator_type.h>
56
57namespace ExaDG
58{
59struct OperatorBaseData
60{
61 OperatorBaseData()
62 : dof_index(0),
63 dof_index_inhomogeneous(dealii::numbers::invalid_unsigned_int),
64 quad_index(0),
65 use_matrix_based_operator_level(false),
66 sparse_matrix_type(SparseMatrixType::Undefined),
67 operator_is_singular(false),
68 use_cell_based_loops(false),
69 implement_block_diagonal_preconditioner_matrix_free(false),
70 solver_block_diagonal(Elementwise::Solver::GMRES),
71 preconditioner_block_diagonal(Elementwise::Preconditioner::InverseMassMatrix),
72 solver_data_block_diagonal(SolverData(1000, 1.e-12, 1.e-2, LinearSolver::Undefined, 1000))
73 {
74 }
75
76 unsigned int dof_index;
77
78 // In addition to the "standard" dof_index above, we need a separate dof index to evaluate
79 // inhomogeneous boundary data correctly. This "inhomogeneous" dof index corresponds to an
80 // AffineConstraints object that only applies periodicity and hanging node constraints.
81 unsigned int dof_index_inhomogeneous;
82
83 unsigned int quad_index;
84
85 // this parameter can be used to use sparse matrices for the vmult() operation. The default
86 // case is to use a matrix-free implementation, i.e. use_matrix_based_operator_level = false.
87 bool use_matrix_based_operator_level;
88
89 SparseMatrixType sparse_matrix_type;
90
91 // Solution of linear systems of equations and preconditioning
92 bool operator_is_singular;
93
94 bool use_cell_based_loops;
95
96 // block Jacobi preconditioner
97 bool implement_block_diagonal_preconditioner_matrix_free;
98
99 // elementwise iterative solution of block Jacobi problems
100 Elementwise::Solver solver_block_diagonal;
101 Elementwise::Preconditioner preconditioner_block_diagonal;
102 SolverData solver_data_block_diagonal;
103};
104
105template<int dim, typename Number, int n_components = 1>
106class OperatorBase : public dealii::EnableObserverPointer
107{
108public:
109 typedef OperatorBase<dim, Number, n_components> This;
110
111 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
112 typedef std::pair<unsigned int, unsigned int> Range;
113 typedef CellIntegrator<dim, n_components, Number> IntegratorCell;
114 typedef FaceIntegrator<dim, n_components, Number> IntegratorFace;
115
116 static unsigned int const dimension = dim;
117 static unsigned int const vectorization_length = dealii::VectorizedArray<Number>::size();
118
119 typedef dealii::LAPACKFullMatrix<Number> LAPACKMatrix;
120
121 typedef dealii::FullMatrix<dealii::TrilinosScalar> FullMatrix_;
122
123 OperatorBase();
124
125 virtual ~OperatorBase()
126 {
127 if(data.sparse_matrix_type == SparseMatrixType::PETSc)
128 {
129#ifdef DEAL_II_WITH_PETSC
130 if(system_matrix_petsc.m() > 0)
131 {
132 PetscErrorCode ierr = VecDestroy(&petsc_vector_dst);
133 AssertThrow(ierr == 0, dealii::ExcPETScError(ierr));
134 ierr = VecDestroy(&petsc_vector_src);
135 AssertThrow(ierr == 0, dealii::ExcPETScError(ierr));
136 }
137#endif
138 }
139 }
140
141 /*
142 * Getters and setters.
143 */
144 void
145 set_time(double const time) const;
146
147 double
148 get_time() const;
149
150 unsigned int
151 get_level() const;
152
153 dealii::AffineConstraints<Number> const &
154 get_affine_constraints() const;
155
156 dealii::MatrixFree<dim, Number> const &
157 get_matrix_free() const;
158
159 unsigned int
160 get_dof_index() const;
161
162 unsigned int
163 get_dof_index_inhomogeneous() const;
164
165 unsigned int
166 get_quad_index() const;
167
168 /*
169 * Returns whether the operator is singular, e.g., the Laplace operator with pure Neumann boundary
170 * conditions is singular.
171 */
172 bool
173 operator_is_singular() const;
174
175 void
176 vmult(VectorType & dst, VectorType const & src) const;
177
178 void
179 vmult_add(VectorType & dst, VectorType const & src) const;
180
181 void
182 vmult_interface_down(VectorType & dst, VectorType const & src) const;
183
184 void
185 vmult_add_interface_up(VectorType & dst, VectorType const & src) const;
186
187 dealii::types::global_dof_index
188 m() const;
189
190 dealii::types::global_dof_index
191 n() const;
192
193 Number
194 el(unsigned int const, unsigned int const) const;
195
196 bool
197 is_empty_locally() const;
198
199 void
200 initialize_dof_vector(VectorType & vector) const;
201
202 /*
203 * For time-dependent problems, the function set_time() needs to be called prior to the present
204 * function.
205 */
206 virtual void
207 set_inhomogeneous_constrained_values(VectorType & solution) const;
208
209 void
210 set_constrained_dofs_to_zero(VectorType & vector) const;
211
212 void
213 calculate_inverse_diagonal(VectorType & diagonal) const;
214
215 /*
216 * Initialize block diagonal preconditioner: initialize everything related to block diagonal
217 * preconditioner.
218 */
219 void
220 initialize_block_diagonal_preconditioner(bool const initialize) const;
221
222 /*
223 * Update block diagonal preconditioner: Recompute block matrices in case of matrix-based
224 * implementation.
225 */
226 void
227 update_block_diagonal_preconditioner() const;
228
233 void
234 apply_inverse_block_diagonal(VectorType & dst, VectorType const & src) const;
235
236 /*
237 * Algebraic multigrid (AMG): sparse matrix (Trilinos) methods
238 */
239#ifdef DEAL_II_WITH_TRILINOS
240 void
241 init_system_matrix(dealii::TrilinosWrappers::SparseMatrix & system_matrix,
242 MPI_Comm const & mpi_comm) const;
243
244 void
245 calculate_system_matrix(dealii::TrilinosWrappers::SparseMatrix & system_matrix) const;
246#endif
247
248 /*
249 * Algebraic multigrid (AMG): sparse matrix (PETSc) methods
250 */
251#ifdef DEAL_II_WITH_PETSC
252 void
253 init_system_matrix(dealii::PETScWrappers::MPI::SparseMatrix & system_matrix,
254 MPI_Comm const & mpi_comm) const;
255
256 void
257 calculate_system_matrix(dealii::PETScWrappers::MPI::SparseMatrix & system_matrix) const;
258#endif
259
260 /*
261 * Provide near null space basis vectors used e.g. in AMG setup. ExaDG assumes a scalar Laplace
262 * operator as default, filling `constant_modes`, which can be overwritten in derived classes if
263 * necessary. `constant_modes_values` may alternatively be used to provide non-trivial modes.
264 */
265 virtual void
266 get_constant_modes(std::vector<std::vector<bool>> & constant_modes,
267 std::vector<std::vector<double>> & constant_modes_values) const;
268
269 /*
270 * Evaluate the homogeneous part of an operator. The homogeneous operator is the operator that is
271 * obtained for homogeneous boundary conditions. This operation is typically applied in linear
272 * iterative solvers (as well as multigrid preconditioners and smoothers). Operations of this type
273 * are called apply_...() and vmult_...() as required by deal.II interfaces.
274 */
275 void
276 apply(VectorType & dst, VectorType const & src) const;
277
278 /*
279 * See function apply() for a description.
280 */
281 void
282 apply_add(VectorType & dst, VectorType const & src) const;
283
284 void
285 assemble_matrix_if_matrix_based() const;
286
287 /*
288 * Matrix-based version of the apply function. This function is used if
289 * use_matrix_based_operator_level = true.
290 */
291 void
292 apply_matrix_based(VectorType & dst, VectorType const & src) const;
293
294 /*
295 * See function apply_matrix_based() for a description.
296 */
297 void
298 apply_matrix_based_add(VectorType & dst, VectorType const & src) const;
299
300 /*
301 * evaluate inhomogeneous parts of operator related to inhomogeneous boundary face integrals.
302 * Operations of this type are called rhs_...() since these functions are called to calculate the
303 * vector forming the right-hand side vector of linear systems of equations. Functions of type
304 * rhs only make sense for linear operators (but they have e.g. no meaning for linearized
305 * operators of nonlinear problems). For this reason, these functions are currently defined
306 * 'virtual' to provide the opportunity to override and assert these functions in derived classes.
307 *
308 * For continuous Galerkin discretizations, this function calls internally the member function
309 * set_inhomogeneous_constrained_values(). Hence, prior to calling this function, one needs to
310 * call set_time() for a correct evaluation in case of time-dependent problems.
311 *
312 * This function sets the dst vector to zero, and afterwards calls rhs_add().
313 */
314 virtual void
315 rhs(VectorType & dst) const;
316
317 /*
318 * See function rhs() for a description.
319 */
320 virtual void
321 rhs_add(VectorType & dst) const;
322
323 /*
324 * Evaluate the operator including homogeneous and inhomogeneous contributions. The typical use
325 * case would be explicit time integration where a splitting into homogeneous and inhomogeneous
326 * contributions is not required. Functions of type evaluate only make sense for linear operators
327 * (but they have e.g. no meaning for linearized operators of nonlinear problems). For this
328 * reason, these functions are currently defined 'virtual' to provide the opportunity to override
329 * and assert these functions in derived classes.
330 *
331 * Unlike the function rhs(), this function does not internally call the function
332 * set_inhomogeneous_constrained_values() prior to evaluation. Hence, one needs to explicitly call
333 * the function set_inhomogeneous_constrained_values() in case of continuous Galerkin
334 * discretizations with inhomogeneous Dirichlet boundary conditions before calling the present
335 * function.
336 */
337 virtual void
338 evaluate(VectorType & dst, VectorType const & src) const;
339
340 /*
341 * See function evaluate() for a description.
342 */
343 virtual void
344 evaluate_add(VectorType & dst, VectorType const & src) const;
345
346 /*
347 * point Jacobi preconditioner (diagonal)
348 */
349 void
350 calculate_diagonal(VectorType & diagonal) const;
351
352 void
353 add_diagonal(VectorType & diagonal) const;
354
355 /*
356 * block Jacobi preconditioner (block-diagonal)
357 */
358 void
359 add_block_diagonal_matrices(std::vector<LAPACKMatrix> & matrices) const;
360
361 void
362 apply_inverse_block_diagonal_matrix_based(VectorType & dst, VectorType const & src) const;
363
364 // matrix-free implementation
365
366 // This function has to initialize everything related to the block diagonal preconditioner when
367 // using the matrix-free variant with elementwise iterative solvers and matrix-free operator
368 // evaluation.
369 void
370 initialize_block_diagonal_preconditioner_matrix_free(bool const initialize) const;
371
372 // updates block diagonal preconditioner when using the matrix-free variant with elementwise
373 // iterative solvers and matrix-free operator evaluation.
374 void
375 update_block_diagonal_preconditioner_matrix_free() const;
376
377 // initializes block diagonal preconditioner for matrix-based variant
378 void
379 initialize_block_diagonal_preconditioner_matrix_based(bool const initialize) const;
380
381 // updates block diagonal preconditioner for matrix-based variant
382 void
383 update_block_diagonal_preconditioner_matrix_based() const;
384
385 void
386 apply_add_block_diagonal_elementwise(unsigned int const cell,
387 dealii::VectorizedArray<Number> * const dst,
388 dealii::VectorizedArray<Number> const * const src,
389 unsigned int const problem_size) const;
390
391 /*
392 * additive Schwarz preconditioner (cellwise block-diagonal)
393 */
394 virtual void
395 compute_factorized_additive_schwarz_matrices() const;
396
397 void
398 apply_inverse_additive_schwarz_matrices(VectorType & dst, VectorType const & src) const;
399
400protected:
401 void
402 reinit(dealii::MatrixFree<dim, Number> const & matrix_free,
403 dealii::AffineConstraints<Number> const & constraints,
404 OperatorBaseData const & data);
405
406 void
407 reinit_cell(IntegratorCell & integrator, unsigned int const cell) const;
408
409 void
410 reinit_face(IntegratorFace & integrator_m,
411 IntegratorFace & integrator_p,
412 unsigned int const face) const;
413
414 void
415 reinit_boundary_face(IntegratorFace & integrator_m, unsigned int const face) const;
416
417 void
418 reinit_face_cell_based(IntegratorFace & integrator_m,
419 IntegratorFace & integrator_p,
420 unsigned int const cell,
421 unsigned int const face,
422 dealii::types::boundary_id const boundary_id) const;
423
424 /*
425 * These methods have to be overwritten by derived classes because these functions are
426 * operator-specific and define how the operator looks like.
427 */
428
429 // standard integration procedure with separate loops for cell and face integrals
430 virtual void
431 do_cell_integral(IntegratorCell & integrator) const;
432
433 virtual void
434 do_face_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const;
435
436 virtual void
437 do_boundary_integral(IntegratorFace & integrator,
438 OperatorType const & operator_type,
439 dealii::types::boundary_id const & boundary_id) const;
440
441 virtual void
442 do_boundary_integral_continuous(IntegratorFace & integrator,
443 OperatorType const & operator_type,
444 dealii::types::boundary_id const & boundary_id) const;
445
446 // The computation of the diagonal and block-diagonal requires face integrals of type
447 // interior (int) and exterior (ext)
448 virtual void
449 do_face_int_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const;
450
451 virtual void
452 do_face_ext_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const;
453
454 // This function is currently only needed due to limitations of deal.II which do
455 // currently not allow to access neighboring data in case of cell-based face loops.
456 // Once this functionality is available, this function should be removed again.
457 // Since only special operators need to evaluate neighboring data, this function
458 // simply redirects to do_face_int_integral() if this function is not overwritten
459 // by a derived class (such as convective terms that require an additional
460 // evaluation of velocity fields for example).
461 virtual void
462 do_face_int_integral_cell_based(IntegratorFace & integrator_m,
463 IntegratorFace & integrator_p) const;
464
465 /*
466 * Matrix-free object.
467 */
469
470 /*
471 * Physical time (required for time-dependent problems).
472 */
473 mutable double time;
474
475 /*
476 * Constraint matrix.
477 */
479
480 mutable dealii::AffineConstraints<double> constraint_double;
481
482 /*
483 * Cell and face integrator flags.
484 */
485 IntegratorFlags integrator_flags;
486
487private:
488 /*
489 * Is the operator used as a multigrid level operator?
490 */
491 bool is_mg;
492
493protected:
494 /*
495 * Is the discretization based on discontinuous Galerkin method?
496 */
497 bool is_dg;
498
499 /*
500 * Block Jacobi preconditioner/smoother: matrix-free version with elementwise iterative solver
501 */
502 typedef Elementwise::OperatorBase<dim, Number, This> ElementwiseOperator;
504 ElementwisePreconditionerBase;
505 typedef Elementwise::
506 IterativeSolver<dim, n_components, Number, ElementwiseOperator, ElementwisePreconditionerBase>
507 ElementwiseSolver;
508
509 mutable std::shared_ptr<ElementwiseOperator> elementwise_operator;
510 mutable std::shared_ptr<ElementwisePreconditionerBase> elementwise_preconditioner;
511 mutable std::shared_ptr<ElementwiseSolver> elementwise_solver;
512
513private:
514 virtual void
515 reinit_cell_derived(IntegratorCell & integrator, unsigned int const cell) const;
516
517 virtual void
518 reinit_face_derived(IntegratorFace & integrator_m,
519 IntegratorFace & integrator_p,
520 unsigned int const face) const;
521
522 virtual void
523 reinit_boundary_face_derived(IntegratorFace & integrator_m, unsigned int const face) const;
524
525 virtual void
526 reinit_face_cell_based_derived(IntegratorFace & integrator_m,
527 IntegratorFace & integrator_p,
528 unsigned int const cell,
529 unsigned int const face,
530 dealii::types::boundary_id const boundary_id) const;
531
532 /*
533 * Helper functions:
534 *
535 * The diagonal, block-diagonal, as well as the system matrix (assembled into a sparse matrix) are
536 * computed columnwise. This means that column i of the block-matrix is computed by evaluating the
537 * operator for a unit vector which takes a value of 1 in row i and is 0 for all other entries.
538 */
539 void
540 create_standard_basis(unsigned int j, IntegratorCell & integrator) const;
541
542 void
543 create_standard_basis(unsigned int j, IntegratorFace & integrator) const;
544
545 void
546 create_standard_basis(unsigned int j,
547 IntegratorFace & integrator_1,
548 IntegratorFace & integrator_2) const;
549
550 /*
551 * This function calculates cell integrals for the full (= homogeneous + inhomogeneous) part,
552 * where inhomogeneous part may occur for continuous Galerkin discretizations due to inhomogeneous
553 * Dirichlet BCs. A prerequisite to call this function is to set inhomogeneous Dirichlet degrees
554 * of freedom appropriately in the DoF-vector src. In case the DoF-vector src is zero apart from
555 * the inhomogeneous Dirichlet degrees of freedom, this function calculates only the inhomogeneous
556 * part of the operator, because the homogeneous operator is zero in this case. For DG
557 * discretizations, this function is equivalent to cell_loop().
558 */
559 void
560 cell_loop_full_operator(dealii::MatrixFree<dim, Number> const & matrix_free,
561 VectorType & dst,
562 VectorType const & src,
563 Range const & range) const;
564
565 /*
566 * This function loops over all cells and calculates cell integrals.
567 */
568 void
569 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
570 VectorType & dst,
571 VectorType const & src,
572 Range const & range) const;
573
574 /*
575 * This function loops over all interior faces and calculates face integrals.
576 */
577 void
578 face_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
579 VectorType & dst,
580 VectorType const & src,
581 Range const & range) const;
582
583 /*
584 * The following functions loop over all boundary faces and calculate boundary face integrals.
585 * Depending on the operator type, we distinguish between boundary face integrals of type
586 * homogeneous, inhomogeneous, and full.
587 */
588
589 // homogeneous operator
590 void
591 boundary_face_loop_hom_operator(dealii::MatrixFree<dim, Number> const & matrix_free,
592 VectorType & dst,
593 VectorType const & src,
594 Range const & range) const;
595
596 // inhomogeneous operator
597 void
598 boundary_face_loop_inhom_operator(dealii::MatrixFree<dim, Number> const & matrix_free,
599 VectorType & dst,
600 VectorType const & src,
601 Range const & range) const;
602
603 // full operator
604 void
605 boundary_face_loop_full_operator(dealii::MatrixFree<dim, Number> const & matrix_free,
606 VectorType & dst,
607 VectorType const & src,
608 Range const & range) const;
609
610 /*
611 * inhomogeneous operator: For the inhomogeneous operator, we only have to calculate boundary face
612 * integrals. The matrix-free implementation, however, does not offer interfaces for boundary face
613 * integrals only. Hence we have to provide empty functions for cell and interior face integrals.
614 */
615 void
616 cell_loop_empty(dealii::MatrixFree<dim, Number> const & matrix_free,
617 VectorType & dst,
618 VectorType const & src,
619 Range const & range) const;
620
621 void
622 face_loop_empty(dealii::MatrixFree<dim, Number> const & matrix_free,
623 VectorType & dst,
624 VectorType const & src,
625 Range const & range) const;
626
627 /*
628 * Calculate matrix and diagonal using `dealii::MatrixFreeTools`.
629 */
630 void
631 cell_compute_matrix(IntegratorCell & integrator) const;
632
633 void
634 face_compute_matrix(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const;
635
636 void
637 boundary_face_compute_matrix(IntegratorFace & integrator_m) const;
638
639 void
640 cell_based_loop_diagonal(dealii::MatrixFree<dim, Number> const & matrix_free,
641 VectorType & dst,
642 VectorType const & src,
643 Range const & range) const;
644
645 /*
646 * Calculate (assemble) block diagonal.
647 */
648 void
649 cell_loop_block_diagonal(dealii::MatrixFree<dim, Number> const & matrix_free,
650 std::vector<LAPACKMatrix> & matrices,
651 std::vector<LAPACKMatrix> const & src,
652 Range const & range) const;
653
654 void
655 face_loop_block_diagonal(dealii::MatrixFree<dim, Number> const & matrix_free,
656 std::vector<LAPACKMatrix> & matrices,
657 std::vector<LAPACKMatrix> const & src,
658 Range const & range) const;
659
660 void
661 boundary_face_loop_block_diagonal(dealii::MatrixFree<dim, Number> const & matrix_free,
662 std::vector<LAPACKMatrix> & matrices,
663 std::vector<LAPACKMatrix> const & src,
664 Range const & range) const;
665
666 // cell-based variant for computation of both cell and face integrals
667 void
668 cell_based_loop_block_diagonal(dealii::MatrixFree<dim, Number> const & matrix_free,
669 std::vector<LAPACKMatrix> & matrices,
670 std::vector<LAPACKMatrix> const & src,
671 Range const & range) const;
672
673 /*
674 * Apply inverse block diagonal:
675 *
676 * we compute dst = B^{-1} * src
677 */
678 void
679 cell_loop_apply_inverse_block_diagonal_matrix_based(
680 dealii::MatrixFree<dim, Number> const & matrix_free,
681 VectorType & dst,
682 VectorType const & src,
683 Range const & range) const;
684
685 /*
686 * Set up sparse matrix internally for templated matrix type (Trilinos or
687 * PETSc matrices)
688 */
689 template<typename SparseMatrix>
690 void
691 internal_init_system_matrix(SparseMatrix & system_matrix,
692 dealii::DynamicSparsityPattern & dsp,
693 MPI_Comm const & mpi_comm) const;
694
695 template<typename SparseMatrix>
696 void
697 internal_calculate_system_matrix(SparseMatrix & system_matrix) const;
698
699 /*
700 * This function sets entries of the DoF-vector corresponding to constraint DoFs to one.
701 */
702 void
703 set_constrained_dofs_to_one(VectorType & vector) const;
704
705 /*
706 * Do we have to evaluate (boundary) face integrals for this operator? For example, operators
707 * such as the mass operator only involve cell integrals.
708 */
709 bool
710 evaluate_face_integrals() const;
711
712 /*
713 * Compute factorized additive Schwarz matrices.
714 */
715 template<typename SparseMatrix>
716 void
717 internal_compute_factorized_additive_schwarz_matrices() const;
718
719 /*
720 * Data structure containing all operator-specific data.
721 */
722 OperatorBaseData data;
723
724 /*
725 * Multigrid level: 0 <= level <= max_level. If the operator is not used as a multigrid
726 * level operator, this variable takes a value of dealii::numbers::invalid_unsigned_int.
727 */
728 unsigned int level;
729
730 /*
731 * Vector of matrices for block-diagonal preconditioners.
732 */
733 mutable std::vector<LAPACKMatrix> matrices;
734
735 /*
736 * Vector with weights for additive Schwarz preconditioner.
737 */
738 mutable VectorType weights;
739
740 unsigned int n_mpi_processes;
741
742 // sparse matrices for matrix-based vmult
743 mutable bool system_matrix_based_been_initialized;
744
745#ifdef DEAL_II_WITH_TRILINOS
746 mutable dealii::TrilinosWrappers::SparseMatrix system_matrix_trilinos;
747#endif
748
749#ifdef DEAL_II_WITH_PETSC
750 mutable dealii::PETScWrappers::MPI::SparseMatrix system_matrix_petsc;
751
752 // PETSc vector objects to avoid re-allocation in every vmult() operation
753 mutable Vec petsc_vector_src;
754 mutable Vec petsc_vector_dst;
755#endif
756};
757
758} // namespace ExaDG
759
760#endif /* EXADG_OPERATORS_OPERATOR_BASE_H_ */
Definition elementwise_operator.h:34
Definition elementwise_preconditioners.h:41
void apply_inverse_block_diagonal(VectorType &dst, VectorType const &src) const
Definition operator_base.cpp:685
Definition lazy_ptr.h:29
Definition driver.cpp:33
Definition integrator_flags.h:31
Definition operator_base.h:60
Definition solver_data.h:92