Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
step-55.h
Go to the documentation of this file.
1
234 *  
235 *   namespace LA
236 *   {
237 *   #if defined(DEAL_II_WITH_PETSC) && !defined(DEAL_II_PETSC_WITH_COMPLEX) && \
238 *   !(defined(DEAL_II_WITH_TRILINOS) && defined(FORCE_USE_OF_TRILINOS))
239 *   using namespace dealii::LinearAlgebraPETSc;
240 *   # define USE_PETSC_LA
241 *   #elif defined(DEAL_II_WITH_TRILINOS)
242 *   using namespace dealii::LinearAlgebraTrilinos;
243 *   #else
244 *   # error DEAL_II_WITH_PETSC or DEAL_II_WITH_TRILINOS required
245 *   #endif
246 *   } // namespace LA
247 *  
248 *   #include <deal.II/lac/vector.h>
249 *   #include <deal.II/lac/full_matrix.h>
250 *   #include <deal.II/lac/solver_cg.h>
251 *   #include <deal.II/lac/solver_gmres.h>
252 *   #include <deal.II/lac/solver_minres.h>
253 *   #include <deal.II/lac/affine_constraints.h>
254 *   #include <deal.II/lac/dynamic_sparsity_pattern.h>
255 *  
256 *   #include <deal.II/lac/petsc_sparse_matrix.h>
257 *   #include <deal.II/lac/petsc_vector.h>
258 *   #include <deal.II/lac/petsc_solver.h>
259 *   #include <deal.II/lac/petsc_precondition.h>
260 *  
261 *   #include <deal.II/grid/grid_generator.h>
262 *   #include <deal.II/grid/manifold_lib.h>
263 *   #include <deal.II/grid/grid_tools.h>
264 *   #include <deal.II/dofs/dof_handler.h>
265 *   #include <deal.II/dofs/dof_renumbering.h>
266 *   #include <deal.II/dofs/dof_tools.h>
267 *   #include <deal.II/fe/fe_values.h>
268 *   #include <deal.II/fe/fe_q.h>
269 *   #include <deal.II/fe/fe_system.h>
270 *   #include <deal.II/numerics/vector_tools.h>
271 *   #include <deal.II/numerics/data_out.h>
272 *   #include <deal.II/numerics/error_estimator.h>
273 *  
274 *   #include <deal.II/base/utilities.h>
275 *   #include <deal.II/base/conditional_ostream.h>
276 *   #include <deal.II/base/index_set.h>
277 *   #include <deal.II/lac/sparsity_tools.h>
278 *   #include <deal.II/distributed/tria.h>
279 *   #include <deal.II/distributed/grid_refinement.h>
280 *  
281 *   #include <cmath>
282 *   #include <fstream>
283 *   #include <iostream>
284 *  
285 *   namespace Step55
286 *   {
287 *   using namespace dealii;
288 *  
289 * @endcode
290 *
291 *
292 * <a name="Linearsolversandpreconditioners"></a>
293 * <h3>Linear solvers and preconditioners</h3>
294 *
295
296 *
297 * We need a few helper classes to represent our solver strategy
298 * described in the introduction.
299 *
300
301 *
302 *
303 * @code
304 *   namespace LinearSolvers
305 *   {
306 * @endcode
307 *
308 * This class exposes the action of applying the inverse of a giving
309 * matrix via the function InverseMatrix::vmult(). Internally, the
310 * inverse is not formed explicitly. Instead, a linear solver with CG
311 * is performed. This class extends the InverseMatrix class in @ref step_22 "step-22"
312 * with an option to specify a preconditioner, and to allow for different
313 * vector types in the vmult function.
314 *
315 * @code
316 *   template <class Matrix, class Preconditioner>
317 *   class InverseMatrix : public Subscriptor
318 *   {
319 *   public:
320 *   InverseMatrix(const Matrix &m, const Preconditioner &preconditioner);
321 *  
322 *   template <typename VectorType>
323 *   void vmult(VectorType &dst, const VectorType &src) const;
324 *  
325 *   private:
327 *   const Preconditioner & preconditioner;
328 *   };
329 *  
330 *  
331 *   template <class Matrix, class Preconditioner>
332 *   InverseMatrix<Matrix, Preconditioner>::InverseMatrix(
333 *   const Matrix & m,
334 *   const Preconditioner &preconditioner)
335 *   : matrix(&m)
336 *   , preconditioner(preconditioner)
337 *   {}
338 *  
339 *  
340 *  
341 *   template <class Matrix, class Preconditioner>
342 *   template <typename VectorType>
343 *   void
344 *   InverseMatrix<Matrix, Preconditioner>::vmult(VectorType & dst,
345 *   const VectorType &src) const
346 *   {
347 *   SolverControl solver_control(src.size(), 1e-8 * src.l2_norm());
348 *   SolverCG<VectorType> cg(solver_control);
349 *   dst = 0;
350 *  
351 *   try
352 *   {
353 *   cg.solve(*matrix, dst, src, preconditioner);
354 *   }
355 *   catch (std::exception &e)
356 *   {
357 *   Assert(false, ExcMessage(e.what()));
358 *   }
359 *   }
360 *  
361 *  
362 * @endcode
363 *
364 * The class A template class for a simple block diagonal preconditioner
365 * for 2x2 matrices.
366 *
367 * @code
368 *   template <class PreconditionerA, class PreconditionerS>
369 *   class BlockDiagonalPreconditioner : public Subscriptor
370 *   {
371 *   public:
372 *   BlockDiagonalPreconditioner(const PreconditionerA &preconditioner_A,
373 *   const PreconditionerS &preconditioner_S);
374 *  
375 *   void vmult(LA::MPI::BlockVector & dst,
376 *   const LA::MPI::BlockVector &src) const;
377 *  
378 *   private:
379 *   const PreconditionerA &preconditioner_A;
380 *   const PreconditionerS &preconditioner_S;
381 *   };
382 *  
383 *   template <class PreconditionerA, class PreconditionerS>
384 *   BlockDiagonalPreconditioner<PreconditionerA, PreconditionerS>::
385 *   BlockDiagonalPreconditioner(const PreconditionerA &preconditioner_A,
386 *   const PreconditionerS &preconditioner_S)
387 *   : preconditioner_A(preconditioner_A)
388 *   , preconditioner_S(preconditioner_S)
389 *   {}
390 *  
391 *  
392 *   template <class PreconditionerA, class PreconditionerS>
393 *   void BlockDiagonalPreconditioner<PreconditionerA, PreconditionerS>::vmult(
394 *   LA::MPI::BlockVector & dst,
395 *   const LA::MPI::BlockVector &src) const
396 *   {
397 *   preconditioner_A.vmult(dst.block(0), src.block(0));
398 *   preconditioner_S.vmult(dst.block(1), src.block(1));
399 *   }
400 *  
401 *   } // namespace LinearSolvers
402 *  
403 * @endcode
404 *
405 *
406 * <a name="Problemsetup"></a>
407 * <h3>Problem setup</h3>
408 *
409
410 *
411 * The following classes represent the right hand side and the exact
412 * solution for the test problem.
413 *
414
415 *
416 *
417 * @code
418 *   template <int dim>
419 *   class RightHandSide : public Function<dim>
420 *   {
421 *   public:
422 *   RightHandSide()
423 *   : Function<dim>(dim + 1)
424 *   {}
425 *  
426 *   virtual void vector_value(const Point<dim> &p,
427 *   Vector<double> & value) const override;
428 *   };
429 *  
430 *  
431 *   template <int dim>
432 *   void RightHandSide<dim>::vector_value(const Point<dim> &p,
433 *   Vector<double> & values) const
434 *   {
435 *   const double R_x = p[0];
436 *   const double R_y = p[1];
437 *  
438 *   constexpr double pi = numbers::PI;
439 *   constexpr double pi2 = numbers::PI * numbers::PI;
440 *  
441 *   values[0] = -1.0L / 2.0L * (-2 * std::sqrt(25.0 + 4 * pi2) + 10.0) *
442 *   std::exp(R_x * (-2 * std::sqrt(25.0 + 4 * pi2) + 10.0)) -
443 *   0.4 * pi2 * std::exp(R_x * (-std::sqrt(25.0 + 4 * pi2) + 5.0)) *
444 *   std::cos(2 * R_y * pi) +
445 *   0.1 * std::pow(-std::sqrt(25.0 + 4 * pi2) + 5.0, 2) *
446 *   std::exp(R_x * (-std::sqrt(25.0 + 4 * pi2) + 5.0)) *
447 *   std::cos(2 * R_y * pi);
448 *   values[1] = 0.2 * pi * (-std::sqrt(25.0 + 4 * pi2) + 5.0) *
449 *   std::exp(R_x * (-std::sqrt(25.0 + 4 * pi2) + 5.0)) *
450 *   std::sin(2 * R_y * pi) -
451 *   0.05 * std::pow(-std::sqrt(25.0 + 4 * pi2) + 5.0, 3) *
452 *   std::exp(R_x * (-std::sqrt(25.0 + 4 * pi2) + 5.0)) *
453 *   std::sin(2 * R_y * pi) / pi;
454 *   values[2] = 0;
455 *   }
456 *  
457 *  
458 *   template <int dim>
459 *   class ExactSolution : public Function<dim>
460 *   {
461 *   public:
462 *   ExactSolution()
463 *   : Function<dim>(dim + 1)
464 *   {}
465 *  
466 *   virtual void vector_value(const Point<dim> &p,
467 *   Vector<double> & values) const override;
468 *   };
469 *  
470 *   template <int dim>
471 *   void ExactSolution<dim>::vector_value(const Point<dim> &p,
472 *   Vector<double> & values) const
473 *   {
474 *   const double R_x = p[0];
475 *   const double R_y = p[1];
476 *  
477 *   constexpr double pi = numbers::PI;
478 *   constexpr double pi2 = numbers::PI * numbers::PI;
479 *  
480 *   values[0] = -std::exp(R_x * (-std::sqrt(25.0 + 4 * pi2) + 5.0)) *
481 *   std::cos(2 * R_y * pi) +
482 *   1;
483 *   values[1] = (1.0L / 2.0L) * (-std::sqrt(25.0 + 4 * pi2) + 5.0) *
484 *   std::exp(R_x * (-std::sqrt(25.0 + 4 * pi2) + 5.0)) *
485 *   std::sin(2 * R_y * pi) / pi;
486 *   values[2] =
487 *   -1.0L / 2.0L * std::exp(R_x * (-2 * std::sqrt(25.0 + 4 * pi2) + 10.0)) -
488 *   2.0 *
489 *   (-6538034.74494422 +
490 *   0.0134758939981709 * std::exp(4 * std::sqrt(25.0 + 4 * pi2))) /
491 *   (-80.0 * std::exp(3 * std::sqrt(25.0 + 4 * pi2)) +
492 *   16.0 * std::sqrt(25.0 + 4 * pi2) *
493 *   std::exp(3 * std::sqrt(25.0 + 4 * pi2))) -
494 *   1634508.68623606 * std::exp(-3.0 * std::sqrt(25.0 + 4 * pi2)) /
495 *   (-10.0 + 2.0 * std::sqrt(25.0 + 4 * pi2)) +
496 *   (-0.00673794699908547 * std::exp(std::sqrt(25.0 + 4 * pi2)) +
497 *   3269017.37247211 * std::exp(-3 * std::sqrt(25.0 + 4 * pi2))) /
498 *   (-8 * std::sqrt(25.0 + 4 * pi2) + 40.0) +
499 *   0.00336897349954273 * std::exp(1.0 * std::sqrt(25.0 + 4 * pi2)) /
500 *   (-10.0 + 2.0 * std::sqrt(25.0 + 4 * pi2));
501 *   }
502 *  
503 *  
504 *  
505 * @endcode
506 *
507 *
508 * <a name="Themainprogram"></a>
509 * <h3>The main program</h3>
510 *
511
512 *
513 * The main class is very similar to @ref step_40 "step-40", except that matrices and
514 * vectors are now block versions, and we store a std::vector<IndexSet>
515 * for owned and relevant DoFs instead of a single IndexSet. We have
516 * exactly two IndexSets, one for all velocity unknowns and one for all
517 * pressure unknowns.
518 *
519 * @code
520 *   template <int dim>
521 *   class StokesProblem
522 *   {
523 *   public:
524 *   StokesProblem(unsigned int velocity_degree);
525 *  
526 *   void run();
527 *  
528 *   private:
529 *   void make_grid();
530 *   void setup_system();
531 *   void assemble_system();
532 *   void solve();
533 *   void refine_grid();
534 *   void output_results(const unsigned int cycle) const;
535 *  
536 *   unsigned int velocity_degree;
537 *   double viscosity;
538 *   MPI_Comm mpi_communicator;
539 *  
540 *   FESystem<dim> fe;
542 *   DoFHandler<dim> dof_handler;
543 *  
544 *   std::vector<IndexSet> owned_partitioning;
545 *   std::vector<IndexSet> relevant_partitioning;
546 *  
547 *   AffineConstraints<double> constraints;
548 *  
549 *   LA::MPI::BlockSparseMatrix system_matrix;
550 *   LA::MPI::BlockSparseMatrix preconditioner_matrix;
551 *   LA::MPI::BlockVector locally_relevant_solution;
552 *   LA::MPI::BlockVector system_rhs;
553 *  
554 *   ConditionalOStream pcout;
555 *   TimerOutput computing_timer;
556 *   };
557 *  
558 *  
559 *  
560 *   template <int dim>
561 *   StokesProblem<dim>::StokesProblem(unsigned int velocity_degree)
562 *   : velocity_degree(velocity_degree)
563 *   , viscosity(0.1)
564 *   , mpi_communicator(MPI_COMM_WORLD)
565 *   , fe(FE_Q<dim>(velocity_degree), dim, FE_Q<dim>(velocity_degree - 1), 1)
566 *   , triangulation(mpi_communicator,
570 *   , dof_handler(triangulation)
571 *   , pcout(std::cout,
572 *   (Utilities::MPI::this_mpi_process(mpi_communicator) == 0))
573 *   , computing_timer(mpi_communicator,
574 *   pcout,
577 *   {}
578 *  
579 *  
580 * @endcode
581 *
582 * The Kovasznay flow is defined on the domain [-0.5, 1.5]^2, which we
583 * create by passing the min and max values to GridGenerator::hyper_cube.
584 *
585 * @code
586 *   template <int dim>
587 *   void StokesProblem<dim>::make_grid()
588 *   {
590 *   triangulation.refine_global(3);
591 *   }
592 *  
593 * @endcode
594 *
595 *
596 * <a name="SystemSetup"></a>
597 * <h3>System Setup</h3>
598 *
599
600 *
601 * The construction of the block matrices and vectors is new compared to
602 * @ref step_40 "step-40" and is different compared to serial codes like @ref step_22 "step-22", because
603 * we need to supply the set of rows that belong to our processor.
604 *
605 * @code
606 *   template <int dim>
607 *   void StokesProblem<dim>::setup_system()
608 *   {
609 *   TimerOutput::Scope t(computing_timer, "setup");
610 *  
611 *   dof_handler.distribute_dofs(fe);
612 *  
613 * @endcode
614 *
615 * Put all dim velocities into block 0 and the pressure into block 1,
616 * then reorder the unknowns by block. Finally count how many unknowns
617 * we have per block.
618 *
619 * @code
620 *   std::vector<unsigned int> stokes_sub_blocks(dim + 1, 0);
621 *   stokes_sub_blocks[dim] = 1;
622 *   DoFRenumbering::component_wise(dof_handler, stokes_sub_blocks);
623 *  
624 *   const std::vector<types::global_dof_index> dofs_per_block =
625 *   DoFTools::count_dofs_per_fe_block(dof_handler, stokes_sub_blocks);
626 *  
627 *   const unsigned int n_u = dofs_per_block[0];
628 *   const unsigned int n_p = dofs_per_block[1];
629 *  
630 *   pcout << " Number of degrees of freedom: " << dof_handler.n_dofs() << " ("
631 *   << n_u << '+' << n_p << ')' << std::endl;
632 *  
633 * @endcode
634 *
635 * We split up the IndexSet for locally owned and locally relevant DoFs
636 * into two IndexSets based on how we want to create the block matrices
637 * and vectors.
638 *
639 * @code
640 *   owned_partitioning.resize(2);
641 *   owned_partitioning[0] = dof_handler.locally_owned_dofs().get_view(0, n_u);
642 *   owned_partitioning[1] =
643 *   dof_handler.locally_owned_dofs().get_view(n_u, n_u + n_p);
644 *  
645 *   const IndexSet locally_relevant_dofs =
647 *   relevant_partitioning.resize(2);
648 *   relevant_partitioning[0] = locally_relevant_dofs.get_view(0, n_u);
649 *   relevant_partitioning[1] = locally_relevant_dofs.get_view(n_u, n_u + n_p);
650 *  
651 * @endcode
652 *
653 * Setting up the constraints for boundary conditions and hanging nodes
654 * is identical to @ref step_40 "step-40". Even though we don't have any hanging nodes
655 * because we only perform global refinement, it is still a good idea
656 * to put this function call in, in case adaptive refinement gets
657 * introduced later.
658 *
659 * @code
660 *   {
661 *   constraints.reinit(locally_relevant_dofs);
662 *  
663 *   const FEValuesExtractors::Vector velocities(0);
664 *   DoFTools::make_hanging_node_constraints(dof_handler, constraints);
665 *   VectorTools::interpolate_boundary_values(dof_handler,
666 *   0,
667 *   ExactSolution<dim>(),
668 *   constraints,
669 *   fe.component_mask(velocities));
670 *   constraints.close();
671 *   }
672 *  
673 * @endcode
674 *
675 * Now we create the system matrix based on a BlockDynamicSparsityPattern.
676 * We know that we won't have coupling between different velocity
677 * components (because we use the laplace and not the deformation tensor)
678 * and no coupling between pressure with its test functions, so we use
679 * a Table to communicate this coupling information to
681 *
682 * @code
683 *   {
684 *   system_matrix.clear();
685 *  
686 *   Table<2, DoFTools::Coupling> coupling(dim + 1, dim + 1);
687 *   for (unsigned int c = 0; c < dim + 1; ++c)
688 *   for (unsigned int d = 0; d < dim + 1; ++d)
689 *   if (c == dim && d == dim)
690 *   coupling[c][d] = DoFTools::none;
691 *   else if (c == dim || d == dim || c == d)
692 *   coupling[c][d] = DoFTools::always;
693 *   else
694 *   coupling[c][d] = DoFTools::none;
695 *  
696 *   BlockDynamicSparsityPattern dsp(relevant_partitioning);
697 *  
699 *   dof_handler, coupling, dsp, constraints, false);
700 *  
702 *   dsp,
703 *   dof_handler.locally_owned_dofs(),
704 *   mpi_communicator,
705 *   locally_relevant_dofs);
706 *  
707 *   system_matrix.reinit(owned_partitioning, dsp, mpi_communicator);
708 *   }
709 *  
710 * @endcode
711 *
712 * The preconditioner matrix has a different coupling (we only fill in
713 * the 1,1 block with the @ref GlossMassMatrix "mass matrix"), otherwise this code is identical
714 * to the construction of the system_matrix above.
715 *
716 * @code
717 *   {
718 *   preconditioner_matrix.clear();
719 *  
720 *   Table<2, DoFTools::Coupling> coupling(dim + 1, dim + 1);
721 *   for (unsigned int c = 0; c < dim + 1; ++c)
722 *   for (unsigned int d = 0; d < dim + 1; ++d)
723 *   if (c == dim && d == dim)
724 *   coupling[c][d] = DoFTools::always;
725 *   else
726 *   coupling[c][d] = DoFTools::none;
727 *  
728 *   BlockDynamicSparsityPattern dsp(relevant_partitioning);
729 *  
731 *   dof_handler, coupling, dsp, constraints, false);
733 *   dsp,
734 *   Utilities::MPI::all_gather(mpi_communicator,
735 *   dof_handler.locally_owned_dofs()),
736 *   mpi_communicator,
737 *   locally_relevant_dofs);
738 *   preconditioner_matrix.reinit(owned_partitioning, dsp, mpi_communicator);
739 *   }
740 *  
741 * @endcode
742 *
743 * Finally, we construct the block vectors with the right sizes. The
744 * function call with two std::vector<IndexSet> will create a ghosted
745 * vector.
746 *
747 * @code
748 *   locally_relevant_solution.reinit(owned_partitioning,
749 *   relevant_partitioning,
750 *   mpi_communicator);
751 *   system_rhs.reinit(owned_partitioning, mpi_communicator);
752 *   }
753 *  
754 *  
755 *  
756 * @endcode
757 *
758 *
759 * <a name="Assembly"></a>
760 * <h3>Assembly</h3>
761 *
762
763 *
764 * This function assembles the system matrix, the preconditioner matrix,
765 * and the right hand side. The code is pretty standard.
766 *
767 * @code
768 *   template <int dim>
769 *   void StokesProblem<dim>::assemble_system()
770 *   {
771 *   TimerOutput::Scope t(computing_timer, "assembly");
772 *  
773 *   system_matrix = 0;
774 *   preconditioner_matrix = 0;
775 *   system_rhs = 0;
776 *  
777 *   const QGauss<dim> quadrature_formula(velocity_degree + 1);
778 *  
779 *   FEValues<dim> fe_values(fe,
780 *   quadrature_formula,
783 *  
784 *   const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
785 *   const unsigned int n_q_points = quadrature_formula.size();
786 *  
787 *   FullMatrix<double> cell_matrix(dofs_per_cell, dofs_per_cell);
788 *   FullMatrix<double> cell_matrix2(dofs_per_cell, dofs_per_cell);
789 *   Vector<double> cell_rhs(dofs_per_cell);
790 *  
791 *   const RightHandSide<dim> right_hand_side;
792 *   std::vector<Vector<double>> rhs_values(n_q_points, Vector<double>(dim + 1));
793 *  
794 *   std::vector<Tensor<2, dim>> grad_phi_u(dofs_per_cell);
795 *   std::vector<double> div_phi_u(dofs_per_cell);
796 *   std::vector<double> phi_p(dofs_per_cell);
797 *  
798 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
799 *   const FEValuesExtractors::Vector velocities(0);
800 *   const FEValuesExtractors::Scalar pressure(dim);
801 *  
802 *   for (const auto &cell : dof_handler.active_cell_iterators())
803 *   if (cell->is_locally_owned())
804 *   {
805 *   cell_matrix = 0;
806 *   cell_matrix2 = 0;
807 *   cell_rhs = 0;
808 *  
809 *   fe_values.reinit(cell);
810 *   right_hand_side.vector_value_list(fe_values.get_quadrature_points(),
811 *   rhs_values);
812 *   for (unsigned int q = 0; q < n_q_points; ++q)
813 *   {
814 *   for (unsigned int k = 0; k < dofs_per_cell; ++k)
815 *   {
816 *   grad_phi_u[k] = fe_values[velocities].gradient(k, q);
817 *   div_phi_u[k] = fe_values[velocities].divergence(k, q);
818 *   phi_p[k] = fe_values[pressure].value(k, q);
819 *   }
820 *  
821 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
822 *   {
823 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
824 *   {
825 *   cell_matrix(i, j) +=
826 *   (viscosity *
827 *   scalar_product(grad_phi_u[i], grad_phi_u[j]) -
828 *   div_phi_u[i] * phi_p[j] - phi_p[i] * div_phi_u[j]) *
829 *   fe_values.JxW(q);
830 *  
831 *   cell_matrix2(i, j) += 1.0 / viscosity * phi_p[i] *
832 *   phi_p[j] * fe_values.JxW(q);
833 *   }
834 *  
835 *   const unsigned int component_i =
836 *   fe.system_to_component_index(i).first;
837 *   cell_rhs(i) += fe_values.shape_value(i, q) *
838 *   rhs_values[q](component_i) * fe_values.JxW(q);
839 *   }
840 *   }
841 *  
842 *  
843 *   cell->get_dof_indices(local_dof_indices);
844 *   constraints.distribute_local_to_global(cell_matrix,
845 *   cell_rhs,
846 *   local_dof_indices,
847 *   system_matrix,
848 *   system_rhs);
849 *  
850 *   constraints.distribute_local_to_global(cell_matrix2,
851 *   local_dof_indices,
852 *   preconditioner_matrix);
853 *   }
854 *  
855 *   system_matrix.compress(VectorOperation::add);
856 *   preconditioner_matrix.compress(VectorOperation::add);
857 *   system_rhs.compress(VectorOperation::add);
858 *   }
859 *  
860 *  
861 *  
862 * @endcode
863 *
864 *
865 * <a name="Solving"></a>
866 * <h3>Solving</h3>
867 *
868
869 *
870 * This function solves the linear system with MINRES with a block diagonal
871 * preconditioner and AMG for the two diagonal blocks as described in the
872 * introduction. The preconditioner applies a v cycle to the 0,0 block
873 * and a CG with the mass matrix for the 1,1 block (the Schur complement).
874 *
875 * @code
876 *   template <int dim>
877 *   void StokesProblem<dim>::solve()
878 *   {
879 *   TimerOutput::Scope t(computing_timer, "solve");
880 *  
881 *   LA::MPI::PreconditionAMG prec_A;
882 *   {
883 *   LA::MPI::PreconditionAMG::AdditionalData data;
884 *  
885 *   #ifdef USE_PETSC_LA
886 *   data.symmetric_operator = true;
887 *   #endif
888 *   prec_A.initialize(system_matrix.block(0, 0), data);
889 *   }
890 *  
891 *   LA::MPI::PreconditionAMG prec_S;
892 *   {
893 *   LA::MPI::PreconditionAMG::AdditionalData data;
894 *  
895 *   #ifdef USE_PETSC_LA
896 *   data.symmetric_operator = true;
897 *   #endif
898 *   prec_S.initialize(preconditioner_matrix.block(1, 1), data);
899 *   }
900 *  
901 * @endcode
902 *
903 * The InverseMatrix is used to solve for the mass matrix:
904 *
905 * @code
906 *   using mp_inverse_t = LinearSolvers::InverseMatrix<LA::MPI::SparseMatrix,
907 *   LA::MPI::PreconditionAMG>;
908 *   const mp_inverse_t mp_inverse(preconditioner_matrix.block(1, 1), prec_S);
909 *  
910 * @endcode
911 *
912 * This constructs the block preconditioner based on the preconditioners
913 * for the individual blocks defined above.
914 *
915 * @code
916 *   const LinearSolvers::BlockDiagonalPreconditioner<LA::MPI::PreconditionAMG,
917 *   mp_inverse_t>
918 *   preconditioner(prec_A, mp_inverse);
919 *  
920 * @endcode
921 *
922 * With that, we can finally set up the linear solver and solve the system:
923 *
924 * @code
925 *   SolverControl solver_control(system_matrix.m(),
926 *   1e-10 * system_rhs.l2_norm());
927 *  
928 *   SolverMinRes<LA::MPI::BlockVector> solver(solver_control);
929 *  
930 *   LA::MPI::BlockVector distributed_solution(owned_partitioning,
931 *   mpi_communicator);
932 *  
933 *   constraints.set_zero(distributed_solution);
934 *  
935 *   solver.solve(system_matrix,
936 *   distributed_solution,
937 *   system_rhs,
938 *   preconditioner);
939 *  
940 *   pcout << " Solved in " << solver_control.last_step() << " iterations."
941 *   << std::endl;
942 *  
943 *   constraints.distribute(distributed_solution);
944 *  
945 * @endcode
946 *
947 * Like in @ref step_56 "step-56", we subtract the mean pressure to allow error
948 * computations against our reference solution, which has a mean value
949 * of zero.
950 *
951 * @code
952 *   locally_relevant_solution = distributed_solution;
953 *   const double mean_pressure =
954 *   VectorTools::compute_mean_value(dof_handler,
955 *   QGauss<dim>(velocity_degree + 2),
956 *   locally_relevant_solution,
957 *   dim);
958 *   distributed_solution.block(1).add(-mean_pressure);
959 *   locally_relevant_solution.block(1) = distributed_solution.block(1);
960 *   }
961 *  
962 *  
963 *  
964 * @endcode
965 *
966 *
967 * <a name="Therest"></a>
968 * <h3>The rest</h3>
969 *
970
971 *
972 * The remainder of the code that deals with mesh refinement, output, and
973 * the main loop is pretty standard.
974 *
975 * @code
976 *   template <int dim>
977 *   void StokesProblem<dim>::refine_grid()
978 *   {
979 *   TimerOutput::Scope t(computing_timer, "refine");
980 *  
981 *   triangulation.refine_global();
982 *   }
983 *  
984 *  
985 *  
986 *   template <int dim>
987 *   void StokesProblem<dim>::output_results(const unsigned int cycle) const
988 *   {
989 *   {
990 *   const ComponentSelectFunction<dim> pressure_mask(dim, dim + 1);
991 *   const ComponentSelectFunction<dim> velocity_mask(std::make_pair(0, dim),
992 *   dim + 1);
993 *  
994 *   Vector<double> cellwise_errors(triangulation.n_active_cells());
995 *   QGauss<dim> quadrature(velocity_degree + 2);
996 *  
997 *   VectorTools::integrate_difference(dof_handler,
998 *   locally_relevant_solution,
999 *   ExactSolution<dim>(),
1000 *   cellwise_errors,
1001 *   quadrature,
1003 *   &velocity_mask);
1004 *  
1005 *   const double error_u_l2 =
1007 *   cellwise_errors,
1009 *  
1010 *   VectorTools::integrate_difference(dof_handler,
1011 *   locally_relevant_solution,
1012 *   ExactSolution<dim>(),
1013 *   cellwise_errors,
1014 *   quadrature,
1016 *   &pressure_mask);
1017 *  
1018 *   const double error_p_l2 =
1020 *   cellwise_errors,
1022 *  
1023 *   pcout << "error: u_0: " << error_u_l2 << " p_0: " << error_p_l2
1024 *   << std::endl;
1025 *   }
1026 *  
1027 *  
1028 *   std::vector<std::string> solution_names(dim, "velocity");
1029 *   solution_names.emplace_back("pressure");
1030 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
1031 *   data_component_interpretation(
1033 *   data_component_interpretation.push_back(
1035 *  
1036 *   DataOut<dim> data_out;
1037 *   data_out.attach_dof_handler(dof_handler);
1038 *   data_out.add_data_vector(locally_relevant_solution,
1039 *   solution_names,
1041 *   data_component_interpretation);
1042 *  
1043 *   LA::MPI::BlockVector interpolated;
1044 *   interpolated.reinit(owned_partitioning, MPI_COMM_WORLD);
1045 *   VectorTools::interpolate(dof_handler, ExactSolution<dim>(), interpolated);
1046 *  
1047 *   LA::MPI::BlockVector interpolated_relevant(owned_partitioning,
1048 *   relevant_partitioning,
1049 *   MPI_COMM_WORLD);
1050 *   interpolated_relevant = interpolated;
1051 *   {
1052 *   std::vector<std::string> solution_names(dim, "ref_u");
1053 *   solution_names.emplace_back("ref_p");
1054 *   data_out.add_data_vector(interpolated_relevant,
1055 *   solution_names,
1057 *   data_component_interpretation);
1058 *   }
1059 *  
1060 *  
1061 *   Vector<float> subdomain(triangulation.n_active_cells());
1062 *   for (unsigned int i = 0; i < subdomain.size(); ++i)
1063 *   subdomain(i) = triangulation.locally_owned_subdomain();
1064 *   data_out.add_data_vector(subdomain, "subdomain");
1065 *  
1066 *   data_out.build_patches();
1067 *  
1068 *   data_out.write_vtu_with_pvtu_record(
1069 *   "./", "solution", cycle, mpi_communicator, 2);
1070 *   }
1071 *  
1072 *  
1073 *  
1074 *   template <int dim>
1075 *   void StokesProblem<dim>::run()
1076 *   {
1077 *   #ifdef USE_PETSC_LA
1078 *   pcout << "Running using PETSc." << std::endl;
1079 *   #else
1080 *   pcout << "Running using Trilinos." << std::endl;
1081 *   #endif
1082 *   const unsigned int n_cycles = 5;
1083 *   for (unsigned int cycle = 0; cycle < n_cycles; ++cycle)
1084 *   {
1085 *   pcout << "Cycle " << cycle << ':' << std::endl;
1086 *  
1087 *   if (cycle == 0)
1088 *   make_grid();
1089 *   else
1090 *   refine_grid();
1091 *  
1092 *   setup_system();
1093 *  
1094 *   assemble_system();
1095 *   solve();
1096 *  
1097 *   if (Utilities::MPI::n_mpi_processes(mpi_communicator) <= 32)
1098 *   {
1099 *   TimerOutput::Scope t(computing_timer, "output");
1100 *   output_results(cycle);
1101 *   }
1102 *  
1103 *   computing_timer.print_summary();
1104 *   computing_timer.reset();
1105 *  
1106 *   pcout << std::endl;
1107 *   }
1108 *   }
1109 *   } // namespace Step55
1110 *  
1111 *  
1112 *  
1113 *   int main(int argc, char *argv[])
1114 *   {
1115 *   try
1116 *   {
1117 *   using namespace dealii;
1118 *   using namespace Step55;
1119 *  
1120 *   Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
1121 *  
1122 *   StokesProblem<2> problem(2);
1123 *   problem.run();
1124 *   }
1125 *   catch (std::exception &exc)
1126 *   {
1127 *   std::cerr << std::endl
1128 *   << std::endl
1129 *   << "----------------------------------------------------"
1130 *   << std::endl;
1131 *   std::cerr << "Exception on processing: " << std::endl
1132 *   << exc.what() << std::endl
1133 *   << "Aborting!" << std::endl
1134 *   << "----------------------------------------------------"
1135 *   << std::endl;
1136 *  
1137 *   return 1;
1138 *   }
1139 *   catch (...)
1140 *   {
1141 *   std::cerr << std::endl
1142 *   << std::endl
1143 *   << "----------------------------------------------------"
1144 *   << std::endl;
1145 *   std::cerr << "Unknown exception!" << std::endl
1146 *   << "Aborting!" << std::endl
1147 *   << "----------------------------------------------------"
1148 *   << std::endl;
1149 *   return 1;
1150 *   }
1151 *  
1152 *   return 0;
1153 *   }
1154 * @endcode
1155<a name="Results"></a><h1>Results</h1>
1156
1157
1158As expected from the discussion above, the number of iterations is independent
1159of the number of processors and only very slightly dependent on @f$h@f$:
1160
1161<table>
1162<tr>
1163 <th colspan="2">PETSc</th>
1164 <th colspan="8">number of processors</th>
1165</tr>
1166<tr>
1167 <th>cycle</th>
1168 <th>dofs</th>
1169 <th>1</th>
1170 <th>2</th>
1171 <th>4</th>
1172 <th>8</th>
1173 <th>16</th>
1174 <th>32</th>
1175 <th>64</th>
1176 <th>128</th>
1177</tr>
1178<tr>
1179 <td>0</td>
1180 <td>659</td>
1181 <td>49</td>
1182 <td>49</td>
1183 <td>49</td>
1184 <td>51</td>
1185 <td>51</td>
1186 <td>51</td>
1187 <td>49</td>
1188 <td>49</td>
1189</tr>
1190<tr>
1191 <td>1</td>
1192 <td>2467</td>
1193 <td>52</td>
1194 <td>52</td>
1195 <td>52</td>
1196 <td>52</td>
1197 <td>52</td>
1198 <td>54</td>
1199 <td>54</td>
1200 <td>53</td>
1201</tr>
1202<tr>
1203 <td>2</td>
1204 <td>9539</td>
1205 <td>56</td>
1206 <td>56</td>
1207 <td>56</td>
1208 <td>54</td>
1209 <td>56</td>
1210 <td>56</td>
1211 <td>54</td>
1212 <td>56</td>
1213</tr>
1214<tr>
1215 <td>3</td>
1216 <td>37507</td>
1217 <td>57</td>
1218 <td>57</td>
1219 <td>57</td>
1220 <td>57</td>
1221 <td>57</td>
1222 <td>56</td>
1223 <td>57</td>
1224 <td>56</td>
1225</tr>
1226<tr>
1227 <td>4</td>
1228 <td>148739</td>
1229 <td>58</td>
1230 <td>59</td>
1231 <td>57</td>
1232 <td>59</td>
1233 <td>57</td>
1234 <td>57</td>
1235 <td>57</td>
1236 <td>57</td>
1237</tr>
1238<tr>
1239 <td>5</td>
1240 <td>592387</td>
1241 <td>60</td>
1242 <td>60</td>
1243 <td>59</td>
1244 <td>59</td>
1245 <td>59</td>
1246 <td>59</td>
1247 <td>59</td>
1248 <td>59</td>
1249</tr>
1250<tr>
1251 <td>6</td>
1252 <td>2364419</td>
1253 <td>62</td>
1254 <td>62</td>
1255 <td>61</td>
1256 <td>61</td>
1257 <td>61</td>
1258 <td>61</td>
1259 <td>61</td>
1260 <td>61</td>
1261</tr>
1262</table>
1263
1264<table>
1265<tr>
1266 <th colspan="2">Trilinos</th>
1267 <th colspan="8">number of processors</th>
1268</tr>
1269<tr>
1270 <th>cycle</th>
1271 <th>dofs</th>
1272 <th>1</th>
1273 <th>2</th>
1274 <th>4</th>
1275 <th>8</th>
1276 <th>16</th>
1277 <th>32</th>
1278 <th>64</th>
1279 <th>128</th>
1280</tr>
1281<tr>
1282 <td>0</td>
1283 <td>659</td>
1284 <td>37</td>
1285 <td>37</td>
1286 <td>37</td>
1287 <td>37</td>
1288 <td>37</td>
1289 <td>37</td>
1290 <td>37</td>
1291 <td>37</td>
1292</tr>
1293<tr>
1294 <td>1</td>
1295 <td>2467</td>
1296 <td>92</td>
1297 <td>89</td>
1298 <td>89</td>
1299 <td>82</td>
1300 <td>86</td>
1301 <td>81</td>
1302 <td>78</td>
1303 <td>78</td>
1304</tr>
1305<tr>
1306 <td>2</td>
1307 <td>9539</td>
1308 <td>102</td>
1309 <td>99</td>
1310 <td>96</td>
1311 <td>95</td>
1312 <td>95</td>
1313 <td>88</td>
1314 <td>83</td>
1315 <td>95</td>
1316</tr>
1317<tr>
1318 <td>3</td>
1319 <td>37507</td>
1320 <td>107</td>
1321 <td>105</td>
1322 <td>104</td>
1323 <td>99</td>
1324 <td>100</td>
1325 <td>96</td>
1326 <td>96</td>
1327 <td>90</td>
1328</tr>
1329<tr>
1330 <td>4</td>
1331 <td>148739</td>
1332 <td>112</td>
1333 <td>112</td>
1334 <td>111</td>
1335 <td>111</td>
1336 <td>127</td>
1337 <td>126</td>
1338 <td>115</td>
1339 <td>117</td>
1340</tr>
1341<tr>
1342 <td>5</td>
1343 <td>592387</td>
1344 <td>116</td>
1345 <td>115</td>
1346 <td>114</td>
1347 <td>112</td>
1348 <td>118</td>
1349 <td>120</td>
1350 <td>131</td>
1351 <td>130</td>
1352</tr>
1353<tr>
1354 <td>6</td>
1355 <td>2364419</td>
1356 <td>130</td>
1357 <td>126</td>
1358 <td>120</td>
1359 <td>120</td>
1360 <td>121</td>
1361 <td>122</td>
1362 <td>121</td>
1363 <td>123</td>
1364</tr>
1365</table>
1366
1367While the PETSc results show a constant number of iterations, the iterations
1368increase when using Trilinos. This is likely because of the different settings
1369used for the AMG preconditioner. For performance reasons we do not allow
1370coarsening below a couple thousand unknowns. As the coarse solver is an exact
1371solve (we are using LU by default), a change in number of levels will
1372influence the quality of a V-cycle. Therefore, a V-cycle is closer to an exact
1373solver for smaller problem sizes.
1374
1375<a name="extensions"></a>
1376<a name="Possibilitiesforextensions"></a><h3>Possibilities for extensions</h3>
1377
1378
1379<a name="InvestigateTrilinositerations"></a><h4>Investigate Trilinos iterations</h4>
1380
1381
1382Play with the smoothers, smoothing steps, and other properties for the
1383Trilinos AMG to achieve an optimal preconditioner.
1384
1385<a name="SolvetheOseenprobleminsteadoftheStokessystem"></a><h4>Solve the Oseen problem instead of the Stokes system</h4>
1386
1387
1388This change requires changing the outer solver to GMRES or BiCGStab, because
1389the system is no longer symmetric.
1390
1391You can prescribe the exact flow solution as @f$b@f$ in the convective term @f$b
1392\cdot \nabla u@f$. This should give the same solution as the original problem,
1393if you set the right hand side to zero.
1394
1395<a name="Adaptiverefinement"></a><h4>Adaptive refinement</h4>
1396
1397
1398So far, this tutorial program refines the mesh globally in each step.
1399Replacing the code in StokesProblem::refine_grid() by something like
1400@code
1401Vector<float> estimated_error_per_cell(triangulation.n_active_cells());
1402
1403FEValuesExtractors::Vector velocities(0);
1404KellyErrorEstimator<dim>::estimate(
1405 dof_handler,
1406 QGauss<dim - 1>(fe.degree + 1),
1407 std::map<types::boundary_id, const Function<dim> *>(),
1408 locally_relevant_solution,
1409 estimated_error_per_cell,
1410 fe.component_mask(velocities));
1411parallel::distributed::GridRefinement::refine_and_coarsen_fixed_number(
1412 triangulation, estimated_error_per_cell, 0.3, 0.0);
1413triangulation.execute_coarsening_and_refinement();
1414@endcode
1415makes it simple to explore adaptive mesh refinement.
1416 *
1417 *
1418<a name="PlainProg"></a>
1419<h1> The plain program</h1>
1420@include "step-55.cc"
1421*/
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
Definition fe_q.h:551
virtual void vector_value(const Point< dim > &p, Vector< RangeNumberType > &values) const
IndexSet get_view(const size_type begin, const size_type end) const
Definition index_set.cc:214
Definition point.h:112
@ wall_times
Definition timer.h:652
__global__ void set(Number *val, const Number s, const size_type N)
#define Assert(cond, exc)
void loop(ITERATOR begin, std_cxx20::type_identity_t< ITERATOR > end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, ASSEMBLER &assembler, const LoopControl &lctrl=LoopControl())
Definition loop.h:439
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
@ update_values
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
std::vector< value_type > split(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)
void component_wise(DoFHandler< dim, spacedim > &dof_handler, const std::vector< unsigned int > &target_component=std::vector< unsigned int >())
IndexSet extract_locally_relevant_dofs(const DoFHandler< dim, spacedim > &dof_handler)
std::vector< types::global_dof_index > count_dofs_per_fe_block(const DoFHandler< dim, spacedim > &dof, const std::vector< unsigned int > &target_block=std::vector< unsigned int >())
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
@ matrix
Contents is actually a matrix.
@ symmetric
Matrix is symmetric.
@ diagonal
Matrix is diagonal.
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.)
Definition advection.h:75
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void distribute_sparsity_pattern(DynamicSparsityPattern &dsp, const IndexSet &locally_owned_rows, const MPI_Comm mpi_comm, const IndexSet &locally_relevant_rows)
void call(const std::function< RT()> &function, internal::return_value< RT > &ret_val)
std::vector< unsigned int > serial(const std::vector< unsigned int > &targets, const std::function< RequestType(const unsigned int)> &create_request, const std::function< AnswerType(const unsigned int, const RequestType &)> &answer_request, const std::function< void(const unsigned int, const AnswerType &)> &process_answer, const MPI_Comm comm)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
Definition mpi.cc:150
std::vector< T > all_gather(const MPI_Comm comm, const T &object_to_send)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
Definition mpi.cc:161
VectorType::value_type compute_mean_value(const hp::MappingCollection< dim, spacedim > &mapping_collection, const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim > &q_collection, const VectorType &v, const unsigned int component)
double compute_global_error(const Triangulation< dim, spacedim > &tria, const InVector &cellwise_error, const NormType &norm, const double exponent=2.)
void integrate_difference(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim, typename InVector::value_type > &exact_solution, OutVector &difference, const Quadrature< dim > &q, const NormType &norm, const Function< spacedim, double > *weight=nullptr, const double exponent=2.)
void interpolate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const ComponentMask &component_mask=ComponentMask())
void run(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
int(&) functions(const void *v1, const void *v2)
static constexpr double PI
Definition numbers.h:259
STL namespace.
::VectorizedArray< Number, width > exp(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
Definition types.h:33
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
const TriangulationDescription::Settings settings
const ::Triangulation< dim, spacedim > & tria