deal.II version 9.7.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
time_dependent_navier_stokes.h
Go to the documentation of this file.
1
203 *  
204 *   #include <deal.II/base/function.h>
205 *   #include <deal.II/base/logstream.h>
206 *   #include <deal.II/base/quadrature_lib.h>
207 *   #include <deal.II/base/quadrature_point_data.h>
208 *   #include <deal.II/base/tensor.h>
209 *   #include <deal.II/base/timer.h>
210 *   #include <deal.II/base/utilities.h>
211 *  
212 *   #include <deal.II/lac/affine_constraints.h>
213 *   #include <deal.II/lac/block_sparse_matrix.h>
214 *   #include <deal.II/lac/block_vector.h>
215 *   #include <deal.II/lac/dynamic_sparsity_pattern.h>
216 *   #include <deal.II/lac/full_matrix.h>
217 *   #include <deal.II/lac/precondition.h>
218 *   #include <deal.II/lac/solver_cg.h>
219 *   #include <deal.II/lac/solver_gmres.h>
220 *   #include <deal.II/lac/sparse_direct.h>
221 *   #include <deal.II/lac/sparsity_tools.h>
222 *  
223 *   #include <deal.II/lac/petsc_block_sparse_matrix.h>
224 *   #include <deal.II/lac/petsc_sparse_matrix.h>
225 *   #include <deal.II/lac/petsc_vector.h>
226 *   #include <deal.II/lac/petsc_precondition.h>
227 *   #include <deal.II/lac/petsc_solver.h>
228 *  
229 *   #include <deal.II/grid/grid_generator.h>
230 *   #include <deal.II/grid/grid_refinement.h>
231 *   #include <deal.II/grid/grid_tools.h>
232 *   #include <deal.II/grid/manifold_lib.h>
233 *   #include <deal.II/grid/tria.h>
234 *   #include <deal.II/grid/tria_accessor.h>
235 *   #include <deal.II/grid/tria_iterator.h>
236 *  
237 *   #include <deal.II/dofs/dof_accessor.h>
238 *   #include <deal.II/dofs/dof_handler.h>
239 *   #include <deal.II/dofs/dof_renumbering.h>
240 *   #include <deal.II/dofs/dof_tools.h>
241 *  
242 *   #include <deal.II/fe/fe_q.h>
243 *   #include <deal.II/fe/fe_system.h>
244 *   #include <deal.II/fe/fe_values.h>
245 *  
246 *   #include <deal.II/numerics/data_out.h>
247 *   #include <deal.II/numerics/error_estimator.h>
248 *   #include <deal.II/numerics/matrix_tools.h>
249 *   #include <deal.II/numerics/vector_tools.h>
250 *  
251 *   #include <deal.II/distributed/grid_refinement.h>
252 *   #include <deal.II/distributed/solution_transfer.h>
253 *   #include <deal.II/distributed/tria.h>
254 *  
255 *   #include <fstream>
256 *   #include <iostream>
257 *   #include <sstream>
258 *  
259 *   namespace fluid
260 *   {
261 *   using namespace dealii;
262 *  
263 * @endcode
264 *
265 *
266 * <a name="time_dependent_navier_stokes.cc-Createthetriangulation"></a>
267 * <h3>Create the triangulation</h3>
268 * The code to create triangulation is copied from
269 * [Martin Kronbichler's
270 * code](https://github.com/kronbichler/adaflo/blob/master/tests/flow_past_cylinder.cc)
271 * with very few modifications.
272 *
273
274 *
275 *
276 * <a name="time_dependent_navier_stokes.cc-Helperfunction"></a>
277 * <h4>Helper function</h4>
278 *
279 * @code
280 *   void create_triangulation_2d(Triangulation<2> &tria,
281 *   bool compute_in_2d = true)
282 *   {
283 *   SphericalManifold<2> boundary(Point<2>(0.5, 0.2));
284 *   Triangulation<2> left, middle, right, tmp, tmp2;
285 *   GridGenerator::subdivided_hyper_rectangle(
286 *   left,
287 *   std::vector<unsigned int>({3U, 4U}),
288 *   Point<2>(),
289 *   Point<2>(0.3, 0.41),
290 *   false);
291 *   GridGenerator::subdivided_hyper_rectangle(
292 *   right,
293 *   std::vector<unsigned int>({18U, 4U}),
294 *   Point<2>(0.7, 0),
295 *   Point<2>(2.5, 0.41),
296 *   false);
297 *  
298 * @endcode
299 *
300 * Create middle part first as a hyper shell.
301 *
302 * @code
303 *   GridGenerator::hyper_shell(middle, Point<2>(0.5, 0.2), 0.05, 0.2, 4, true);
304 *   middle.reset_all_manifolds();
305 *   for (Triangulation<2>::cell_iterator cell = middle.begin();
306 *   cell != middle.end();
307 *   ++cell)
308 *   for (unsigned int f = 0; f < GeometryInfo<2>::faces_per_cell; ++f)
309 *   {
310 *   bool is_inner_rim = true;
311 *   for (unsigned int v = 0; v < GeometryInfo<2>::vertices_per_face; ++v)
312 *   {
313 *   Point<2> &vertex = cell->face(f)->vertex(v);
314 *   if (std::abs(vertex.distance(Point<2>(0.5, 0.2)) - 0.05) > 1e-10)
315 *   {
316 *   is_inner_rim = false;
317 *   break;
318 *   }
319 *   }
320 *   if (is_inner_rim)
321 *   cell->face(f)->set_manifold_id(1);
322 *   }
323 *   middle.set_manifold(1, boundary);
324 *   middle.refine_global(1);
325 *  
326 * @endcode
327 *
328 * Then move the vertices to the points where we want them to be to create a
329 * slightly asymmetric cube with a hole:
330 *
331 * @code
332 *   for (Triangulation<2>::cell_iterator cell = middle.begin();
333 *   cell != middle.end();
334 *   ++cell)
335 *   for (unsigned int v = 0; v < GeometryInfo<2>::vertices_per_cell; ++v)
336 *   {
337 *   Point<2> &vertex = cell->vertex(v);
338 *   if (std::abs(vertex[0] - 0.7) < 1e-10 &&
339 *   std::abs(vertex[1] - 0.2) < 1e-10)
340 *   vertex = Point<2>(0.7, 0.205);
341 *   else if (std::abs(vertex[0] - 0.6) < 1e-10 &&
342 *   std::abs(vertex[1] - 0.3) < 1e-10)
343 *   vertex = Point<2>(0.7, 0.41);
344 *   else if (std::abs(vertex[0] - 0.6) < 1e-10 &&
345 *   std::abs(vertex[1] - 0.1) < 1e-10)
346 *   vertex = Point<2>(0.7, 0);
347 *   else if (std::abs(vertex[0] - 0.5) < 1e-10 &&
348 *   std::abs(vertex[1] - 0.4) < 1e-10)
349 *   vertex = Point<2>(0.5, 0.41);
350 *   else if (std::abs(vertex[0] - 0.5) < 1e-10 &&
351 *   std::abs(vertex[1] - 0.0) < 1e-10)
352 *   vertex = Point<2>(0.5, 0.0);
353 *   else if (std::abs(vertex[0] - 0.4) < 1e-10 &&
354 *   std::abs(vertex[1] - 0.3) < 1e-10)
355 *   vertex = Point<2>(0.3, 0.41);
356 *   else if (std::abs(vertex[0] - 0.4) < 1e-10 &&
357 *   std::abs(vertex[1] - 0.1) < 1e-10)
358 *   vertex = Point<2>(0.3, 0);
359 *   else if (std::abs(vertex[0] - 0.3) < 1e-10 &&
360 *   std::abs(vertex[1] - 0.2) < 1e-10)
361 *   vertex = Point<2>(0.3, 0.205);
362 *   else if (std::abs(vertex[0] - 0.56379) < 1e-4 &&
363 *   std::abs(vertex[1] - 0.13621) < 1e-4)
364 *   vertex = Point<2>(0.59, 0.11);
365 *   else if (std::abs(vertex[0] - 0.56379) < 1e-4 &&
366 *   std::abs(vertex[1] - 0.26379) < 1e-4)
367 *   vertex = Point<2>(0.59, 0.29);
368 *   else if (std::abs(vertex[0] - 0.43621) < 1e-4 &&
369 *   std::abs(vertex[1] - 0.13621) < 1e-4)
370 *   vertex = Point<2>(0.41, 0.11);
371 *   else if (std::abs(vertex[0] - 0.43621) < 1e-4 &&
372 *   std::abs(vertex[1] - 0.26379) < 1e-4)
373 *   vertex = Point<2>(0.41, 0.29);
374 *   }
375 *  
376 * @endcode
377 *
378 * Refine once to create the same level of refinement as in the
379 * neighboring domains:
380 *
381 * @code
382 *   middle.refine_global(1);
383 *  
384 * @endcode
385 *
386 * Must copy the triangulation because we cannot merge triangulations with
387 * refinement:
388 *
389 * @code
390 *   GridGenerator::flatten_triangulation(middle, tmp2);
391 *  
392 * @endcode
393 *
394 * Left domain is required in 3d only.
395 *
396 * @code
397 *   if (compute_in_2d)
398 *   {
399 *   GridGenerator::merge_triangulations(tmp2, right, tria);
400 *   }
401 *   else
402 *   {
403 *   GridGenerator::merge_triangulations(left, tmp2, tmp);
404 *   GridGenerator::merge_triangulations(tmp, right, tria);
405 *   }
406 *   }
407 *  
408 * @endcode
409 *
410 *
411 * <a name="time_dependent_navier_stokes.cc-2Dflowaroundcylindertriangulation"></a>
412 * <h4>2D flow around cylinder triangulation</h4>
413 *
414 * @code
415 *   void create_triangulation(Triangulation<2> &tria)
416 *   {
417 *   create_triangulation_2d(tria);
418 * @endcode
419 *
420 * Set the left boundary (inflow) to 0, the right boundary (outflow) to 1,
421 * upper to 2, lower to 3 and the cylindrical surface to 4.
422 *
423 * @code
424 *   for (Triangulation<2>::active_cell_iterator cell = tria.begin();
425 *   cell != tria.end();
426 *   ++cell)
427 *   {
428 *   for (unsigned int f = 0; f < GeometryInfo<2>::faces_per_cell; ++f)
429 *   {
430 *   if (cell->face(f)->at_boundary())
431 *   {
432 *   if (std::abs(cell->face(f)->center()[0] - 2.5) < 1e-12)
433 *   {
434 *   cell->face(f)->set_all_boundary_ids(1);
435 *   }
436 *   else if (std::abs(cell->face(f)->center()[0] - 0.3) < 1e-12)
437 *   {
438 *   cell->face(f)->set_all_boundary_ids(0);
439 *   }
440 *   else if (std::abs(cell->face(f)->center()[1] - 0.41) < 1e-12)
441 *   {
442 *   cell->face(f)->set_all_boundary_ids(3);
443 *   }
444 *   else if (std::abs(cell->face(f)->center()[1]) < 1e-12)
445 *   {
446 *   cell->face(f)->set_all_boundary_ids(2);
447 *   }
448 *   else
449 *   {
450 *   cell->face(f)->set_all_boundary_ids(4);
451 *   }
452 *   }
453 *   }
454 *   }
455 *   }
456 *  
457 * @endcode
458 *
459 *
460 * <a name="time_dependent_navier_stokes.cc-3Dflowaroundcylindertriangulation"></a>
461 * <h4>3D flow around cylinder triangulation</h4>
462 *
463 * @code
464 *   void create_triangulation(Triangulation<3> &tria)
465 *   {
466 *   Triangulation<2> tria_2d;
467 *   create_triangulation_2d(tria_2d, false);
468 *   GridGenerator::extrude_triangulation(tria_2d, 5, 0.41, tria);
469 * @endcode
470 *
471 * Set the ids of the boundaries in x direction to 0 and 1; y direction to 2 and 3;
472 * z direction to 4 and 5; the cylindrical surface 6.
473 *
474 * @code
475 *   for (Triangulation<3>::active_cell_iterator cell = tria.begin();
476 *   cell != tria.end();
477 *   ++cell)
478 *   {
479 *   for (unsigned int f = 0; f < GeometryInfo<3>::faces_per_cell; ++f)
480 *   {
481 *   if (cell->face(f)->at_boundary())
482 *   {
483 *   if (std::abs(cell->face(f)->center()[0] - 2.5) < 1e-12)
484 *   {
485 *   cell->face(f)->set_all_boundary_ids(1);
486 *   }
487 *   else if (std::abs(cell->face(f)->center()[0]) < 1e-12)
488 *   {
489 *   cell->face(f)->set_all_boundary_ids(0);
490 *   }
491 *   else if (std::abs(cell->face(f)->center()[1] - 0.41) < 1e-12)
492 *   {
493 *   cell->face(f)->set_all_boundary_ids(3);
494 *   }
495 *   else if (std::abs(cell->face(f)->center()[1]) < 1e-12)
496 *   {
497 *   cell->face(f)->set_all_boundary_ids(2);
498 *   }
499 *   else if (std::abs(cell->face(f)->center()[2] - 0.41) < 1e-12)
500 *   {
501 *   cell->face(f)->set_all_boundary_ids(5);
502 *   }
503 *   else if (std::abs(cell->face(f)->center()[2]) < 1e-12)
504 *   {
505 *   cell->face(f)->set_all_boundary_ids(4);
506 *   }
507 *   else
508 *   {
509 *   cell->face(f)->set_all_boundary_ids(6);
510 *   }
511 *   }
512 *   }
513 *   }
514 *   }
515 *  
516 * @endcode
517 *
518 *
519 * <a name="time_dependent_navier_stokes.cc-Timestepping"></a>
520 * <h3>Time stepping</h3>
521 * This class is pretty much self-explanatory.
522 *
523 * @code
524 *   class Time
525 *   {
526 *   public:
527 *   Time(const double time_end,
528 *   const double delta_t,
529 *   const double output_interval,
530 *   const double refinement_interval)
531 *   : timestep(0),
532 *   time_current(0.0),
533 *   time_end(time_end),
534 *   delta_t(delta_t),
535 *   output_interval(output_interval),
536 *   refinement_interval(refinement_interval)
537 *   {
538 *   }
539 *   double current() const { return time_current; }
540 *   double end() const { return time_end; }
541 *   double get_delta_t() const { return delta_t; }
542 *   unsigned int get_timestep() const { return timestep; }
543 *   bool time_to_output() const;
544 *   bool time_to_refine() const;
545 *   void increment();
546 *  
547 *   private:
548 *   unsigned int timestep;
549 *   double time_current;
550 *   const double time_end;
551 *   const double delta_t;
552 *   const double output_interval;
553 *   const double refinement_interval;
554 *   };
555 *  
556 *   bool Time::time_to_output() const
557 *   {
558 *   unsigned int delta = static_cast<unsigned int>(output_interval / delta_t);
559 *   return (timestep >= delta && timestep % delta == 0);
560 *   }
561 *  
562 *   bool Time::time_to_refine() const
563 *   {
564 *   unsigned int delta = static_cast<unsigned int>(refinement_interval / delta_t);
565 *   return (timestep >= delta && timestep % delta == 0);
566 *   }
567 *  
568 *   void Time::increment()
569 *   {
570 *   time_current += delta_t;
571 *   ++timestep;
572 *   }
573 *  
574 * @endcode
575 *
576 *
577 * <a name="time_dependent_navier_stokes.cc-Boundaryvalues"></a>
578 * <h3>Boundary values</h3>
579 * Dirichlet boundary conditions for the velocity inlet and walls.
580 *
581 * @code
582 *   template <int dim>
583 *   class BoundaryValues : public Function<dim>
584 *   {
585 *   public:
586 *   BoundaryValues() : Function<dim>(dim + 1) {}
587 *   virtual double value(const Point<dim> &p,
588 *   const unsigned int component) const override;
589 *  
590 *   virtual void vector_value(const Point<dim> &p,
591 *   Vector<double> &values) const override;
592 *   };
593 *  
594 *   template <int dim>
595 *   double BoundaryValues<dim>::value(const Point<dim> &p,
596 *   const unsigned int component) const
597 *   {
598 *   Assert(component < this->n_components,
599 *   ExcIndexRange(component, 0, this->n_components));
600 *   double left_boundary = (dim == 2 ? 0.3 : 0.0);
601 *   if (component == 0 && std::abs(p[0] - left_boundary) < 1e-10)
602 *   {
603 * @endcode
604 *
605 * For a parabolic velocity profile, @f$U_\mathrm{avg} = 2/3
606 * U_\mathrm{max}@f$
607 * in 2D, and @f$U_\mathrm{avg} = 4/9 U_\mathrm{max}@f$ in 3D.
608 * If @f$\nu = 0.001@f$, @f$D = 0.1@f$, then @f$Re = 100 U_\mathrm{avg}@f$.
609 *
610 * @code
611 *   double Uavg = 1.0;
612 *   double Umax = (dim == 2 ? 3 * Uavg / 2 : 9 * Uavg / 4);
613 *   double value = 4 * Umax * p[1] * (0.41 - p[1]) / (0.41 * 0.41);
614 *   if (dim == 3)
615 *   {
616 *   value *= 4 * p[2] * (0.41 - p[2]) / (0.41 * 0.41);
617 *   }
618 *   return value;
619 *   }
620 *   return 0;
621 *   }
622 *  
623 *   template <int dim>
624 *   void BoundaryValues<dim>::vector_value(const Point<dim> &p,
625 *   Vector<double> &values) const
626 *   {
627 *   for (unsigned int c = 0; c < this->n_components; ++c)
628 *   values(c) = BoundaryValues<dim>::value(p, c);
629 *   }
630 *  
631 * @endcode
632 *
633 *
634 * <a name="time_dependent_navier_stokes.cc-Blockpreconditioner"></a>
635 * <h3>Block preconditioner</h3>
636 *
637
638 *
639 * The block Schur preconditioner can be written as the product of three
640 * matrices:
641 * @f$
642 * P^{-1} = \begin{pmatrix} \tilde{A}^{-1} & 0\\ 0 & I\end{pmatrix}
643 * \begin{pmatrix} I & -B^T\\ 0 & I\end{pmatrix}
644 * \begin{pmatrix} I & 0\\ 0 & \tilde{S}^{-1}\end{pmatrix}
645 * @f$
646 * @f$\tilde{A}@f$ is symmetric since the convection term is eliminated from the
647 * LHS.
648 * @f$\tilde{S}^{-1}@f$ is the inverse of the Schur complement of @f$\tilde{A}@f$,
649 * which consists of a reaction term, a diffusion term, a Grad-Div term
650 * and a convection term.
651 * In practice, the convection contribution is ignored, namely
652 * @f$\tilde{S}^{-1} = -(\nu + \gamma)M_p^{-1} -
653 * \frac{1}{\Delta{t}}{[B(diag(M_u))^{-1}B^T]}^{-1}@f$
654 * where @f$M_p@f$ is the pressure mass, and
655 * @f${[B(diag(M_u))^{-1}B^T]}@f$ is an approximation to the Schur complement of
656 * (velocity) mass matrix @f$BM_u^{-1}B^T@f$.
657 *
658
659 *
660 * Same as the tutorials, we define a vmult operation for the block
661 * preconditioner
662 * instead of write it as a matrix. It can be seen from the above definition,
663 * the result of the vmult operation of the block preconditioner can be
664 * obtained
665 * from the results of the vmult operations of @f$M_u^{-1}@f$, @f$M_p^{-1}@f$,
666 * @f$\tilde{A}^{-1}@f$, which can be transformed into solving three symmetric
667 * linear
668 * systems.
669 *
670 * @code
671 *   class BlockSchurPreconditioner : public Subscriptor
672 *   {
673 *   public:
674 *   BlockSchurPreconditioner(
675 *   TimerOutput &timer,
676 *   double gamma,
677 *   double viscosity,
678 *   double dt,
679 *   const std::vector<IndexSet> &owned_partitioning,
680 *   const PETScWrappers::MPI::BlockSparseMatrix &system,
681 *   const PETScWrappers::MPI::BlockSparseMatrix &mass,
682 *   PETScWrappers::MPI::BlockSparseMatrix &schur);
683 *  
684 *   void vmult(PETScWrappers::MPI::BlockVector &dst,
685 *   const PETScWrappers::MPI::BlockVector &src) const;
686 *  
687 *   private:
688 *   TimerOutput &timer;
689 *   const double gamma;
690 *   const double viscosity;
691 *   const double dt;
692 *  
693 *   const SmartPointer<const PETScWrappers::MPI::BlockSparseMatrix>
694 *   system_matrix;
695 *   const SmartPointer<const PETScWrappers::MPI::BlockSparseMatrix> mass_matrix;
696 * @endcode
697 *
698 * As discussed, @f${[B(diag(M_u))^{-1}B^T]}@f$ and its inverse
699 * need to be computed.
700 * We can either explicitly compute it out as a matrix, or define
701 * it as a class with a vmult operation.
702 * The second approach saves some computation to construct the matrix,
703 * but leads to slow convergence in CG solver because it is impossible
704 * to apply a preconditioner. We go with the first route.
705 *
706 * @code
707 *   const SmartPointer<PETScWrappers::MPI::BlockSparseMatrix> mass_schur;
708 *   };
709 *  
710 * @endcode
711 *
712 *
713 * <a name="time_dependent_navier_stokes.cc-BlockSchurPreconditionerBlockSchurPreconditioner"></a>
714 * <h4>BlockSchurPreconditioner::BlockSchurPreconditioner</h4>
715 *
716
717 *
718 * Input parameters and system matrix, mass matrix as well as the mass schur
719 * matrix are needed in the preconditioner. In addition, we pass the
720 * partitioning information into this class because we need to create some
721 * temporary block vectors inside.
722 *
723 * @code
724 *   BlockSchurPreconditioner::BlockSchurPreconditioner(
725 *   TimerOutput &timer,
726 *   double gamma,
727 *   double viscosity,
728 *   double dt,
729 *   const std::vector<IndexSet> &owned_partitioning,
730 *   const PETScWrappers::MPI::BlockSparseMatrix &system,
731 *   const PETScWrappers::MPI::BlockSparseMatrix &mass,
732 *   PETScWrappers::MPI::BlockSparseMatrix &schur)
733 *   : timer(timer),
734 *   gamma(gamma),
735 *   viscosity(viscosity),
736 *   dt(dt),
737 *   system_matrix(&system),
738 *   mass_matrix(&mass),
739 *   mass_schur(&schur)
740 *   {
741 *   TimerOutput::Scope timer_section(timer, "CG for Sm");
742 * @endcode
743 *
744 * The schur complemete of mass matrix is actually being computed here.
745 *
746 * @code
747 *   PETScWrappers::MPI::BlockVector tmp1, tmp2;
748 *   tmp1.reinit(owned_partitioning, mass_matrix->get_mpi_communicator());
749 *   tmp2.reinit(owned_partitioning, mass_matrix->get_mpi_communicator());
750 *   tmp1 = 1;
751 *   tmp2 = 0;
752 * @endcode
753 *
754 * Jacobi preconditioner of matrix A is by definition @f${diag(A)}^{-1}@f$,
755 * this is exactly what we want to compute.
756 *
757 * @code
758 *   PETScWrappers::PreconditionJacobi jacobi(mass_matrix->block(0, 0));
759 *   jacobi.vmult(tmp2.block(0), tmp1.block(0));
760 *   system_matrix->block(1, 0).mmult(
761 *   mass_schur->block(1, 1), system_matrix->block(0, 1), tmp2.block(0));
762 *   }
763 *  
764 * @endcode
765 *
766 *
767 * <a name="time_dependent_navier_stokes.cc-BlockSchurPreconditionervmult"></a>
768 * <h4>BlockSchurPreconditioner::vmult</h4>
769 *
770
771 *
772 * The vmult operation strictly follows the definition of
773 * BlockSchurPreconditioner
774 * introduced above. Conceptually it computes @f$u = P^{-1}v@f$.
775 *
776 * @code
777 *   void BlockSchurPreconditioner::vmult(
778 *   PETScWrappers::MPI::BlockVector &dst,
779 *   const PETScWrappers::MPI::BlockVector &src) const
780 *   {
781 * @endcode
782 *
783 * Temporary vectors
784 *
785 * @code
786 *   PETScWrappers::MPI::Vector utmp(src.block(0));
787 *   PETScWrappers::MPI::Vector tmp(src.block(1));
788 *   tmp = 0;
789 * @endcode
790 *
791 * This block computes @f$u_1 = \tilde{S}^{-1} v_1@f$,
792 * where CG solvers are used for @f$M_p^{-1}@f$ and @f$S_m^{-1}@f$.
793 *
794 * @code
795 *   {
796 *   TimerOutput::Scope timer_section(timer, "CG for Mp");
797 *   SolverControl mp_control(src.block(1).size(),
798 *   1e-6 * src.block(1).l2_norm());
799 *   PETScWrappers::SolverCG cg_mp(mp_control);
800 * @endcode
801 *
802 * @f$-(\nu + \gamma)M_p^{-1}v_1@f$
803 *
804 * @code
805 *   PETScWrappers::PreconditionBlockJacobi Mp_preconditioner;
806 *   Mp_preconditioner.initialize(mass_matrix->block(1, 1));
807 *   cg_mp.solve(
808 *   mass_matrix->block(1, 1), tmp, src.block(1), Mp_preconditioner);
809 *   tmp *= -(viscosity + gamma);
810 *   }
811 * @endcode
812 *
813 * @f$-\frac{1}{dt}S_m^{-1}v_1@f$
814 *
815 * @code
816 *   {
817 *   TimerOutput::Scope timer_section(timer, "CG for Sm");
818 *   SolverControl sm_control(src.block(1).size(),
819 *   1e-6 * src.block(1).l2_norm());
820 *   PETScWrappers::SolverCG cg_sm(sm_control);
821 * @endcode
822 *
823 * PreconditionBlockJacobi works find on Sm if we do not refine the mesh.
824 * Because after refine_mesh is called, zero entries will be created on
825 * the diagonal (not sure why), which prevents PreconditionBlockJacobi
826 * from being used.
827 *
828 * @code
829 *   PETScWrappers::PreconditionNone Sm_preconditioner;
830 *   Sm_preconditioner.initialize(mass_schur->block(1, 1));
831 *   cg_sm.solve(
832 *   mass_schur->block(1, 1), dst.block(1), src.block(1), Sm_preconditioner);
833 *   dst.block(1) *= -1 / dt;
834 *   }
835 * @endcode
836 *
837 * Adding up these two, we get @f$\tilde{S}^{-1}v_1@f$.
838 *
839 * @code
840 *   dst.block(1) += tmp;
841 * @endcode
842 *
843 * Compute @f$v_0 - B^T\tilde{S}^{-1}v_1@f$ based on @f$u_1@f$.
844 *
845 * @code
846 *   system_matrix->block(0, 1).vmult(utmp, dst.block(1));
847 *   utmp *= -1.0;
848 *   utmp += src.block(0);
849 * @endcode
850 *
851 * Finally, compute the product of @f$\tilde{A}^{-1}@f$ and utmp
852 * using another CG solver.
853 *
854 * @code
855 *   {
856 *   TimerOutput::Scope timer_section(timer, "CG for A");
857 *   SolverControl a_control(src.block(0).size(),
858 *   1e-6 * src.block(0).l2_norm());
859 *   PETScWrappers::SolverCG cg_a(a_control);
860 * @endcode
861 *
862 * We do not use any preconditioner for this block, which is of course
863 * slow,
864 * only because the performance of the only two preconditioners available
865 * PreconditionBlockJacobi and PreconditionBoomerAMG are even worse than
866 * none.
867 *
868 * @code
869 *   PETScWrappers::PreconditionNone A_preconditioner;
870 *   A_preconditioner.initialize(system_matrix->block(0, 0));
871 *   cg_a.solve(
872 *   system_matrix->block(0, 0), dst.block(0), utmp, A_preconditioner);
873 *   }
874 *   }
875 *  
876 * @endcode
877 *
878 *
879 * <a name="time_dependent_navier_stokes.cc-TheincompressibleNavierStokessolver"></a>
880 * <h3>The incompressible Navier-Stokes solver</h3>
881 *
882
883 *
884 * Parallel incompressible Navier Stokes equation solver using
885 * implicit-explicit time scheme.
886 * This program is built upon dealii tutorials @ref step_57 "step-57", @ref step_40 "step-40", @ref step_22 "step-22",
887 * and @ref step_20 "step-20".
888 * The system equation is written in the incremental form, and we treat
889 * the convection term explicitly. Therefore the system equation is linear
890 * and symmetric, which does not need to be solved with Newton's iteration.
891 * The system is further stabilized and preconditioned with Grad-Div method,
892 * where GMRES solver is used as the outer solver.
893 *
894 * @code
895 *   template <int dim>
896 *   class InsIMEX
897 *   {
898 *   public:
900 *   void run();
901 *   ~InsIMEX() { timer.print_summary(); }
902 *  
903 *   private:
904 *   void setup_dofs();
905 *   void make_constraints();
906 *   void initialize_system();
907 *   void assemble(bool use_nonzero_constraints, bool assemble_system);
908 *   std::pair<unsigned int, double> solve(bool use_nonzero_constraints,
909 *   bool assemble_system);
910 *   void refine_mesh(const unsigned int, const unsigned int);
911 *   void output_results(const unsigned int) const;
912 *   double viscosity;
913 *   double gamma;
914 *   const unsigned int degree;
915 *   std::vector<types::global_dof_index> dofs_per_block;
916 *  
918 *   FESystem<dim> fe;
919 *   DoFHandler<dim> dof_handler;
920 *   QGauss<dim> volume_quad_formula;
921 *   QGauss<dim - 1> face_quad_formula;
922 *  
923 *   AffineConstraints<double> zero_constraints;
924 *   AffineConstraints<double> nonzero_constraints;
925 *  
926 *   BlockSparsityPattern sparsity_pattern;
927 * @endcode
928 *
929 * System matrix to be solved
930 *
931 * @code
932 *   PETScWrappers::MPI::BlockSparseMatrix system_matrix;
933 * @endcode
934 *
935 * Mass matrix is a block matrix which includes both velocity
936 * mass matrix and pressure mass matrix.
937 *
938 * @code
940 * @endcode
941 *
942 * The schur complement of mass matrix is not a block matrix.
943 * However, because we want to reuse the partition we created
944 * for the system matrix, it is defined as a block matrix
945 * where only one block is actually used.
946 *
947 * @code
949 * @endcode
950 *
951 * The latest known solution.
952 *
953 * @code
954 *   PETScWrappers::MPI::BlockVector present_solution;
955 * @endcode
956 *
957 * The increment at a certain time step.
958 *
959 * @code
960 *   PETScWrappers::MPI::BlockVector solution_increment;
961 * @endcode
962 *
963 * System RHS
964 *
965 * @code
966 *   PETScWrappers::MPI::BlockVector system_rhs;
967 *  
968 *   MPI_Comm mpi_communicator;
969 *  
970 *   ConditionalOStream pcout;
971 *  
972 * @endcode
973 *
974 * The IndexSets of owned velocity and pressure respectively.
975 *
976 * @code
977 *   std::vector<IndexSet> owned_partitioning;
978 *  
979 * @endcode
980 *
981 * The IndexSets of relevant velocity and pressure respectively.
982 *
983 * @code
984 *   std::vector<IndexSet> relevant_partitioning;
985 *  
986 * @endcode
987 *
988 * The IndexSet of all relevant dofs.
989 *
990 * @code
991 *   IndexSet locally_relevant_dofs;
992 *  
993 * @endcode
994 *
995 * The BlockSchurPreconditioner for the entire system.
996 *
997 * @code
998 *   std::shared_ptr<BlockSchurPreconditioner> preconditioner;
999 *  
1000 *   Time time;
1001 *   mutable TimerOutput timer;
1002 *   };
1003 *  
1004 * @endcode
1005 *
1006 *
1007 * <a name="time_dependent_navier_stokes.cc-InsIMEXInsIMEX"></a>
1008 * <h4>InsIMEX::InsIMEX</h4>
1009 *
1010 * @code
1011 *   template <int dim>
1012 *   InsIMEX<dim>::InsIMEX(parallel::distributed::Triangulation<dim> &tria)
1013 *   : viscosity(0.001),
1014 *   gamma(0.1),
1015 *   degree(1),
1016 *   triangulation(tria),
1017 *   fe(FE_Q<dim>(degree + 1), dim, FE_Q<dim>(degree), 1),
1018 *   dof_handler(triangulation),
1019 *   volume_quad_formula(degree + 2),
1020 *   face_quad_formula(degree + 2),
1021 *   mpi_communicator(MPI_COMM_WORLD),
1022 *   pcout(std::cout, Utilities::MPI::this_mpi_process(mpi_communicator) == 0),
1023 *   time(1e0, 1e-3, 1e-2, 1e-2),
1024 *   timer(
1025 *   mpi_communicator, pcout, TimerOutput::never, TimerOutput::wall_times)
1026 *   {
1027 *   }
1028 *  
1029 * @endcode
1030 *
1031 *
1032 * <a name="time_dependent_navier_stokes.cc-InsIMEXsetup_dofs"></a>
1033 * <h4>InsIMEX::setup_dofs</h4>
1034 *
1035 * @code
1036 *   template <int dim>
1037 *   void InsIMEX<dim>::setup_dofs()
1038 *   {
1039 * @endcode
1040 *
1041 * The first step is to associate DoFs with a given mesh.
1042 *
1043 * @code
1044 *   dof_handler.distribute_dofs(fe);
1045 * @endcode
1046 *
1047 * We renumber the components to have all velocity DoFs come before
1048 * the pressure DoFs to be able to split the solution vector in two blocks
1049 * which are separately accessed in the block preconditioner.
1050 *
1051 * @code
1052 *   DoFRenumbering::Cuthill_McKee(dof_handler);
1053 *   std::vector<unsigned int> block_component(dim + 1, 0);
1054 *   block_component[dim] = 1;
1055 *   DoFRenumbering::component_wise(dof_handler, block_component);
1056 *   dofs_per_block = DoFTools::count_dofs_per_fe_block(dof_handler, block_component);
1057 * @endcode
1058 *
1059 * Partitioning.
1060 *
1061 * @code
1062 *   unsigned int dof_u = dofs_per_block[0];
1063 *   unsigned int dof_p = dofs_per_block[1];
1064 *   owned_partitioning.resize(2);
1065 *   owned_partitioning[0] = dof_handler.locally_owned_dofs().get_view(0, dof_u);
1066 *   owned_partitioning[1] =
1067 *   dof_handler.locally_owned_dofs().get_view(dof_u, dof_u + dof_p);
1068 *   locally_relevant_dofs = DoFTools::extract_locally_relevant_dofs(dof_handler);
1069 *   relevant_partitioning.resize(2);
1070 *   relevant_partitioning[0] = locally_relevant_dofs.get_view(0, dof_u);
1071 *   relevant_partitioning[1] =
1072 *   locally_relevant_dofs.get_view(dof_u, dof_u + dof_p);
1073 *   pcout << " Number of active fluid cells: "
1074 *   << triangulation.n_global_active_cells() << std::endl
1075 *   << " Number of degrees of freedom: " << dof_handler.n_dofs() << " ("
1076 *   << dof_u << '+' << dof_p << ')' << std::endl;
1077 *   }
1078 *  
1079 * @endcode
1080 *
1081 *
1082 * <a name="time_dependent_navier_stokes.cc-InsIMEXmake_constraints"></a>
1083 * <h4>InsIMEX::make_constraints</h4>
1084 *
1085 * @code
1086 *   template <int dim>
1087 *   void InsIMEX<dim>::make_constraints()
1088 *   {
1089 * @endcode
1090 *
1091 * Because the equation is written in incremental form, two constraints
1092 * are needed: nonzero constraint and zero constraint.
1093 *
1094 * @code
1095 *   nonzero_constraints.clear();
1096 *   zero_constraints.clear();
1097 *   nonzero_constraints.reinit(locally_relevant_dofs);
1098 *   zero_constraints.reinit(locally_relevant_dofs);
1099 *   DoFTools::make_hanging_node_constraints(dof_handler, nonzero_constraints);
1100 *   DoFTools::make_hanging_node_constraints(dof_handler, zero_constraints);
1101 *  
1102 * @endcode
1103 *
1104 * Apply Dirichlet boundary conditions on all boundaries except for the
1105 * outlet.
1106 *
1107 * @code
1108 *   std::vector<unsigned int> dirichlet_bc_ids;
1109 *   if (dim == 2)
1110 *   dirichlet_bc_ids = std::vector<unsigned int>{0, 2, 3, 4};
1111 *   else
1112 *   dirichlet_bc_ids = std::vector<unsigned int>{0, 2, 3, 4, 5, 6};
1113 *  
1114 *   const FEValuesExtractors::Vector velocities(0);
1115 *   for (auto id : dirichlet_bc_ids)
1116 *   {
1118 *   id,
1119 *   BoundaryValues<dim>(),
1120 *   nonzero_constraints,
1121 *   fe.component_mask(velocities));
1123 *   dof_handler,
1124 *   id,
1125 *   Functions::ZeroFunction<dim>(dim + 1),
1126 *   zero_constraints,
1127 *   fe.component_mask(velocities));
1128 *   }
1129 *   nonzero_constraints.close();
1130 *   zero_constraints.close();
1131 *   }
1132 *  
1133 * @endcode
1134 *
1135 *
1136 * <a name="time_dependent_navier_stokes.cc-InsIMEXinitialize_system"></a>
1137 * <h4>InsIMEX::initialize_system</h4>
1138 *
1139 * @code
1140 *   template <int dim>
1141 *   void InsIMEX<dim>::initialize_system()
1142 *   {
1143 *   preconditioner.reset();
1144 *   system_matrix.clear();
1145 *   mass_matrix.clear();
1146 *   mass_schur.clear();
1147 *  
1148 *   BlockDynamicSparsityPattern dsp(dofs_per_block, dofs_per_block);
1149 *   DoFTools::make_sparsity_pattern(dof_handler, dsp, nonzero_constraints);
1150 *   sparsity_pattern.copy_from(dsp);
1152 *   dsp,
1153 *   dof_handler.locally_owned_dofs(),
1154 *   mpi_communicator,
1155 *   locally_relevant_dofs);
1156 *  
1157 *   system_matrix.reinit(owned_partitioning, dsp, mpi_communicator);
1158 *   mass_matrix.reinit(owned_partitioning, dsp, mpi_communicator);
1159 *  
1160 * @endcode
1161 *
1162 * Only the @f$(1, 1)@f$ block in the mass schur matrix is used.
1163 * Compute the sparsity pattern for mass schur in advance.
1164 * The only nonzero block has the same sparsity pattern as @f$BB^T@f$.
1165 *
1166 * @code
1167 *   BlockDynamicSparsityPattern schur_dsp(dofs_per_block, dofs_per_block);
1168 *   schur_dsp.block(1, 1).compute_mmult_pattern(sparsity_pattern.block(1, 0),
1169 *   sparsity_pattern.block(0, 1));
1170 *   mass_schur.reinit(owned_partitioning, schur_dsp, mpi_communicator);
1171 *  
1172 * @endcode
1173 *
1174 * present_solution is ghosted because it is used in the
1175 * output and mesh refinement functions.
1176 *
1177 * @code
1178 *   present_solution.reinit(
1179 *   owned_partitioning, relevant_partitioning, mpi_communicator);
1180 * @endcode
1181 *
1182 * solution_increment is non-ghosted because the linear solver needs
1183 * a completely distributed vector.
1184 *
1185 * @code
1186 *   solution_increment.reinit(owned_partitioning, mpi_communicator);
1187 * @endcode
1188 *
1189 * system_rhs is non-ghosted because it is only used in the linear
1190 * solver and residual evaluation.
1191 *
1192 * @code
1193 *   system_rhs.reinit(owned_partitioning, mpi_communicator);
1194 *   }
1195 *  
1196 * @endcode
1197 *
1198 *
1199 * <a name="time_dependent_navier_stokes.cc-InsIMEXassemble"></a>
1200 * <h4>InsIMEX::assemble</h4>
1201 *
1202
1203 *
1204 * Assemble the system matrix, mass matrix, and the RHS.
1205 * It can be used to assemble the entire system or only the RHS.
1206 * An additional option is added to determine whether nonzero
1207 * constraints or zero constraints should be used.
1208 * Note that we only need to assemble the LHS for twice: once with the nonzero
1209 * constraint
1210 * and once for zero constraint. But we must assemble the RHS at every time
1211 * step.
1212 *
1213 * @code
1214 *   template <int dim>
1215 *   void InsIMEX<dim>::assemble(bool use_nonzero_constraints,
1216 *   bool assemble_system)
1217 *   {
1218 *   TimerOutput::Scope timer_section(timer, "Assemble system");
1219 *  
1220 *   if (assemble_system)
1221 *   {
1222 *   system_matrix = 0;
1223 *   mass_matrix = 0;
1224 *   }
1225 *   system_rhs = 0;
1226 *  
1227 *   FEValues<dim> fe_values(fe,
1228 *   volume_quad_formula,
1231 *   FEFaceValues<dim> fe_face_values(fe,
1232 *   face_quad_formula,
1235 *   update_JxW_values);
1236 *  
1237 *   const unsigned int dofs_per_cell = fe.dofs_per_cell;
1238 *   const unsigned int n_q_points = volume_quad_formula.size();
1239 *  
1240 *   const FEValuesExtractors::Vector velocities(0);
1241 *   const FEValuesExtractors::Scalar pressure(dim);
1242 *  
1243 *   FullMatrix<double> local_matrix(dofs_per_cell, dofs_per_cell);
1244 *   FullMatrix<double> local_mass_matrix(dofs_per_cell, dofs_per_cell);
1245 *   Vector<double> local_rhs(dofs_per_cell);
1246 *  
1247 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1248 *  
1249 *   std::vector<Tensor<1, dim>> current_velocity_values(n_q_points);
1250 *   std::vector<Tensor<2, dim>> current_velocity_gradients(n_q_points);
1251 *   std::vector<double> current_velocity_divergences(n_q_points);
1252 *   std::vector<double> current_pressure_values(n_q_points);
1253 *  
1254 *   std::vector<double> div_phi_u(dofs_per_cell);
1255 *   std::vector<Tensor<1, dim>> phi_u(dofs_per_cell);
1256 *   std::vector<Tensor<2, dim>> grad_phi_u(dofs_per_cell);
1257 *   std::vector<double> phi_p(dofs_per_cell);
1258 *  
1259 *   for (auto cell = dof_handler.begin_active(); cell != dof_handler.end();
1260 *   ++cell)
1261 *   {
1262 *   if (cell->is_locally_owned())
1263 *   {
1264 *   fe_values.reinit(cell);
1265 *  
1266 *   if (assemble_system)
1267 *   {
1268 *   local_matrix = 0;
1269 *   local_mass_matrix = 0;
1270 *   }
1271 *   local_rhs = 0;
1272 *  
1273 *   fe_values[velocities].get_function_values(present_solution,
1274 *   current_velocity_values);
1275 *  
1276 *   fe_values[velocities].get_function_gradients(
1277 *   present_solution, current_velocity_gradients);
1278 *  
1279 *   fe_values[velocities].get_function_divergences(
1280 *   present_solution, current_velocity_divergences);
1281 *  
1282 *   fe_values[pressure].get_function_values(present_solution,
1283 *   current_pressure_values);
1284 *  
1285 * @endcode
1286 *
1287 * Assemble the system matrix and mass matrix simultaneouly.
1288 * The mass matrix only uses the @f$(0, 0)@f$ and @f$(1, 1)@f$ blocks.
1289 *
1290 * @code
1291 *   for (unsigned int q = 0; q < n_q_points; ++q)
1292 *   {
1293 *   for (unsigned int k = 0; k < dofs_per_cell; ++k)
1294 *   {
1295 *   div_phi_u[k] = fe_values[velocities].divergence(k, q);
1296 *   grad_phi_u[k] = fe_values[velocities].gradient(k, q);
1297 *   phi_u[k] = fe_values[velocities].value(k, q);
1298 *   phi_p[k] = fe_values[pressure].value(k, q);
1299 *   }
1300 *  
1301 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1302 *   {
1303 *   if (assemble_system)
1304 *   {
1305 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
1306 *   {
1307 *   local_matrix(i, j) +=
1308 *   (viscosity *
1309 *   scalar_product(grad_phi_u[j], grad_phi_u[i]) -
1310 *   div_phi_u[i] * phi_p[j] -
1311 *   phi_p[i] * div_phi_u[j] +
1312 *   gamma * div_phi_u[j] * div_phi_u[i] +
1313 *   phi_u[i] * phi_u[j] / time.get_delta_t()) *
1314 *   fe_values.JxW(q);
1315 *   local_mass_matrix(i, j) +=
1316 *   (phi_u[i] * phi_u[j] + phi_p[i] * phi_p[j]) *
1317 *   fe_values.JxW(q);
1318 *   }
1319 *   }
1320 *   local_rhs(i) -=
1321 *   (viscosity * scalar_product(current_velocity_gradients[q],
1322 *   grad_phi_u[i]) -
1323 *   current_velocity_divergences[q] * phi_p[i] -
1324 *   current_pressure_values[q] * div_phi_u[i] +
1325 *   gamma * current_velocity_divergences[q] * div_phi_u[i] +
1326 *   current_velocity_gradients[q] *
1327 *   current_velocity_values[q] * phi_u[i]) *
1328 *   fe_values.JxW(q);
1329 *   }
1330 *   }
1331 *  
1332 *   cell->get_dof_indices(local_dof_indices);
1333 *  
1334 *   const AffineConstraints<double> &constraints_used =
1335 *   use_nonzero_constraints ? nonzero_constraints : zero_constraints;
1336 *   if (assemble_system)
1337 *   {
1338 *   constraints_used.distribute_local_to_global(local_matrix,
1339 *   local_rhs,
1340 *   local_dof_indices,
1341 *   system_matrix,
1342 *   system_rhs);
1343 *   constraints_used.distribute_local_to_global(
1344 *   local_mass_matrix, local_dof_indices, mass_matrix);
1345 *   }
1346 *   else
1347 *   {
1348 *   constraints_used.distribute_local_to_global(
1349 *   local_rhs, local_dof_indices, system_rhs);
1350 *   }
1351 *   }
1352 *   }
1353 *  
1354 *   if (assemble_system)
1355 *   {
1356 *   system_matrix.compress(VectorOperation::add);
1357 *   mass_matrix.compress(VectorOperation::add);
1358 *   }
1359 *   system_rhs.compress(VectorOperation::add);
1360 *   }
1361 *  
1362 * @endcode
1363 *
1364 *
1365 * <a name="time_dependent_navier_stokes.cc-InsIMEXsolve"></a>
1366 * <h4>InsIMEX::solve</h4>
1367 * Solve the linear system using FGMRES solver with block preconditioner.
1368 * After solving the linear system, the same AffineConstraints object as used
1369 * in assembly must be used again, to set the constrained value.
1370 * The second argument is used to determine whether the block
1371 * preconditioner should be reset or not.
1372 *
1373 * @code
1374 *   template <int dim>
1375 *   std::pair<unsigned int, double>
1376 *   InsIMEX<dim>::solve(bool use_nonzero_constraints, bool assemble_system)
1377 *   {
1378 *   if (assemble_system)
1379 *   {
1380 *   preconditioner.reset(new BlockSchurPreconditioner(timer,
1381 *   gamma,
1382 *   viscosity,
1383 *   time.get_delta_t(),
1384 *   owned_partitioning,
1385 *   system_matrix,
1386 *   mass_matrix,
1387 *   mass_schur));
1388 *   }
1389 *  
1390 *   SolverControl solver_control(
1391 *   system_matrix.m(), 1e-8 * system_rhs.l2_norm(), true);
1392 * @endcode
1393 *
1394 * Because PETScWrappers::SolverGMRES only accepts preconditioner
1395 * derived from PETScWrappers::PreconditionBase,
1396 * we use dealii SolverFGMRES.
1397 *
1398 * @code
1400 *   SolverFGMRES<PETScWrappers::MPI::BlockVector> gmres(solver_control,
1401 *   vector_memory);
1402 *  
1403 * @endcode
1404 *
1405 * The solution vector must be non-ghosted
1406 *
1407 * @code
1408 *   gmres.solve(system_matrix, solution_increment, system_rhs, *preconditioner);
1409 *  
1410 *   const AffineConstraints<double> &constraints_used =
1411 *   use_nonzero_constraints ? nonzero_constraints : zero_constraints;
1412 *   constraints_used.distribute(solution_increment);
1413 *  
1414 *   return {solver_control.last_step(), solver_control.last_value()};
1415 *   }
1416 *  
1417 * @endcode
1418 *
1419 *
1420 * <a name="time_dependent_navier_stokes.cc-InsIMEXrun"></a>
1421 * <h4>InsIMEX::run</h4>
1422 *
1423 * @code
1424 *   template <int dim>
1425 *   void InsIMEX<dim>::run()
1426 *   {
1427 *   pcout << "Running with PETSc on "
1428 *   << Utilities::MPI::n_mpi_processes(mpi_communicator)
1429 *   << " MPI rank(s)..." << std::endl;
1430 *  
1431 *   triangulation.refine_global(0);
1432 *   setup_dofs();
1433 *   make_constraints();
1434 *   initialize_system();
1435 *  
1436 * @endcode
1437 *
1438 * Time loop.
1439 *
1440 * @code
1441 *   bool refined = false;
1442 *   while (time.end() - time.current() > 1e-12)
1443 *   {
1444 *   if (time.get_timestep() == 0)
1445 *   {
1446 *   output_results(0);
1447 *   }
1448 *   time.increment();
1449 *   std::cout.precision(6);
1450 *   std::cout.width(12);
1451 *   pcout << std::string(96, '*') << std::endl
1452 *   << "Time step = " << time.get_timestep()
1453 *   << ", at t = " << std::scientific << time.current() << std::endl;
1454 * @endcode
1455 *
1456 * Resetting
1457 *
1458 * @code
1459 *   solution_increment = 0;
1460 * @endcode
1461 *
1462 * Only use nonzero constraints at the very first time step
1463 *
1464 * @code
1465 *   bool apply_nonzero_constraints = (time.get_timestep() == 1);
1466 * @endcode
1467 *
1468 * We have to assemble the LHS for the initial two time steps:
1469 * once using nonzero_constraints, once using zero_constraints,
1470 * as well as the steps immediately after mesh refinement.
1471 *
1472 * @code
1473 *   bool assemble_system = (time.get_timestep() < 3 || refined);
1474 *   refined = false;
1475 *   assemble(apply_nonzero_constraints, assemble_system);
1476 *   auto state = solve(apply_nonzero_constraints, assemble_system);
1477 * @endcode
1478 *
1479 * Note we have to use a non-ghosted vector to do the addition.
1480 *
1481 * @code
1483 *   tmp.reinit(owned_partitioning, mpi_communicator);
1484 *   tmp = present_solution;
1485 *   tmp += solution_increment;
1486 *   present_solution = tmp;
1487 *   pcout << std::scientific << std::left << " GMRES_ITR = " << std::setw(3)
1488 *   << state.first << " GMRES_RES = " << state.second << std::endl;
1489 * @endcode
1490 *
1491 * Output
1492 *
1493 * @code
1494 *   if (time.time_to_output())
1495 *   {
1496 *   output_results(time.get_timestep());
1497 *   }
1498 *   if (time.time_to_refine())
1499 *   {
1500 *   refine_mesh(0, 4);
1501 *   refined = true;
1502 *   }
1503 *   }
1504 *   }
1505 *  
1506 * @endcode
1507 *
1508 *
1509 * <a name="time_dependent_navier_stokes.cc-InsIMEXoutput_result"></a>
1510 * <h4>InsIMEX::output_result</h4>
1511 *
1512
1513 *
1514 *
1515 * @code
1516 *   template <int dim>
1517 *   void InsIMEX<dim>::output_results(const unsigned int output_index) const
1518 *   {
1519 *   TimerOutput::Scope timer_section(timer, "Output results");
1520 *   pcout << "Writing results..." << std::endl;
1521 *   std::vector<std::string> solution_names(dim, "velocity");
1522 *   solution_names.push_back("pressure");
1523 *  
1524 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
1525 *   data_component_interpretation(
1527 *   data_component_interpretation.push_back(
1529 *   DataOut<dim> data_out;
1530 *   data_out.attach_dof_handler(dof_handler);
1531 * @endcode
1532 *
1533 * vector to be output must be ghosted
1534 *
1535 * @code
1536 *   data_out.add_data_vector(present_solution,
1537 *   solution_names,
1539 *   data_component_interpretation);
1540 *  
1541 * @endcode
1542 *
1543 * Partition
1544 *
1545 * @code
1546 *   Vector<float> subdomain(triangulation.n_active_cells());
1547 *   for (unsigned int i = 0; i < subdomain.size(); ++i)
1548 *   {
1549 *   subdomain(i) = triangulation.locally_owned_subdomain();
1550 *   }
1551 *   data_out.add_data_vector(subdomain, "subdomain");
1552 *  
1553 *   data_out.build_patches(degree + 1);
1554 *  
1555 *   std::string basename =
1556 *   "navierstokes" + Utilities::int_to_string(output_index, 6) + "-";
1557 *  
1558 *   std::string filename =
1559 *   basename +
1560 *   Utilities::int_to_string(triangulation.locally_owned_subdomain(), 4) +
1561 *   ".vtu";
1562 *  
1563 *   std::ofstream output(filename);
1564 *   data_out.write_vtu(output);
1565 *  
1566 *   static std::vector<std::pair<double, std::string>> times_and_names;
1567 *   if (Utilities::MPI::this_mpi_process(mpi_communicator) == 0)
1568 *   {
1569 *   for (unsigned int i = 0;
1570 *   i < Utilities::MPI::n_mpi_processes(mpi_communicator);
1571 *   ++i)
1572 *   {
1573 *   times_and_names.push_back(
1574 *   {time.current(),
1575 *   basename + Utilities::int_to_string(i, 4) + ".vtu"});
1576 *   }
1577 *   std::ofstream pvd_output("navierstokes.pvd");
1578 *   DataOutBase::write_pvd_record(pvd_output, times_and_names);
1579 *   }
1580 *   }
1581 *  
1582 * @endcode
1583 *
1584 *
1585 * <a name="time_dependent_navier_stokes.cc-InsIMEXrefine_mesh"></a>
1586 * <h4>InsIMEX::refine_mesh</h4>
1587 *
1588
1589 *
1590 *
1591 * @code
1592 *   template <int dim>
1593 *   void InsIMEX<dim>::refine_mesh(const unsigned int min_grid_level,
1594 *   const unsigned int max_grid_level)
1595 *   {
1596 *   TimerOutput::Scope timer_section(timer, "Refine mesh");
1597 *   pcout << "Refining mesh..." << std::endl;
1598 *  
1599 *   Vector<float> estimated_error_per_cell(triangulation.n_active_cells());
1600 *   const FEValuesExtractors::Vector velocity(0);
1601 *   KellyErrorEstimator<dim>::estimate(dof_handler,
1602 *   face_quad_formula,
1603 *   {},
1604 *   present_solution,
1605 *   estimated_error_per_cell,
1606 *   fe.component_mask(velocity));
1608 *   triangulation, estimated_error_per_cell, 0.6, 0.4);
1609 *   if (triangulation.n_levels() > max_grid_level)
1610 *   {
1611 *   for (auto cell = triangulation.begin_active(max_grid_level);
1612 *   cell != triangulation.end();
1613 *   ++cell)
1614 *   {
1615 *   cell->clear_refine_flag();
1616 *   }
1617 *   }
1618 *   for (auto cell = triangulation.begin_active(min_grid_level);
1619 *   cell != triangulation.end_active(min_grid_level);
1620 *   ++cell)
1621 *   {
1622 *   cell->clear_coarsen_flag();
1623 *   }
1624 *  
1625 * @endcode
1626 *
1627 * Prepare to transfer
1628 *
1629 * @code
1632 *   trans(dof_handler);
1633 *  
1634 *   triangulation.prepare_coarsening_and_refinement();
1635 *  
1636 *   trans.prepare_for_coarsening_and_refinement(present_solution);
1637 *  
1638 * @endcode
1639 *
1640 * Refine the mesh
1641 *
1642 * @code
1643 *   triangulation.execute_coarsening_and_refinement();
1644 *  
1645 * @endcode
1646 *
1647 * Reinitialize the system
1648 *
1649 * @code
1650 *   setup_dofs();
1651 *   make_constraints();
1652 *   initialize_system();
1653 *  
1654 * @endcode
1655 *
1656 * Transfer solution
1657 * Need a non-ghosted vector for interpolation
1658 *
1659 * @code
1660 *   PETScWrappers::MPI::BlockVector tmp(solution_increment);
1661 *   tmp = 0;
1662 *   trans.interpolate(tmp);
1663 *   present_solution = tmp;
1664 *   }
1665 *   }
1666 *  
1667 * @endcode
1668 *
1669 *
1670 * <a name="time_dependent_navier_stokes.cc-mainfunction"></a>
1671 * <h3>main function</h3>
1672 *
1673
1674 *
1675 *
1676 * @code
1677 *   int main(int argc, char *argv[])
1678 *   {
1679 *   try
1680 *   {
1681 *   using namespace dealii;
1682 *   using namespace fluid;
1683 *  
1684 *   Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
1685 *   parallel::distributed::Triangulation<2> tria(MPI_COMM_WORLD);
1686 *   create_triangulation(tria);
1687 *   InsIMEX<2> flow(tria);
1688 *   flow.run();
1689 *   }
1690 *   catch (std::exception &exc)
1691 *   {
1692 *   std::cerr << std::endl
1693 *   << std::endl
1694 *   << "----------------------------------------------------"
1695 *   << std::endl;
1696 *   std::cerr << "Exception on processing: " << std::endl
1697 *   << exc.what() << std::endl
1698 *   << "Aborting!" << std::endl
1699 *   << "----------------------------------------------------"
1700 *   << std::endl;
1701 *   return 1;
1702 *   }
1703 *   catch (...)
1704 *   {
1705 *   std::cerr << std::endl
1706 *   << std::endl
1707 *   << "----------------------------------------------------"
1708 *   << std::endl;
1709 *   std::cerr << "Unknown exception!" << std::endl
1710 *   << "Aborting!" << std::endl
1711 *   << "----------------------------------------------------"
1712 *   << std::endl;
1713 *   return 1;
1714 *   }
1715 *   return 0;
1716 *   }
1717 * @endcode
1718
1719
1720*/
void distribute(VectorType &vec) const
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
Definition fe_q.h:554
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ReadVector< Number > &solution, Vector< float > &error, const ComponentMask &component_mask={}, const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
@ wall_times
Definition timer.h:651
void loop(IteratorType begin, std_cxx20::type_identity_t< IteratorType > end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, AssemblerType &assembler, const LoopControl &lctrl=LoopControl())
Definition loop.h:564
const bool IsBlockVector< VectorType >::value
void make_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ 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)
const Event initial
Definition event.cc:68
void write_pvd_record(std::ostream &out, const std::vector< std::pair< double, std::string > > &times_and_names)
void component_wise(DoFHandler< dim, spacedim > &dof_handler, const std::vector< unsigned int > &target_component=std::vector< unsigned int >())
void Cuthill_McKee(DoFHandler< dim, spacedim > &dof_handler, const bool reversed_numbering=false, const bool use_constraints=false, const std::vector< types::global_dof_index > &starting_indices=std::vector< types::global_dof_index >())
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 create_triangulation(Triangulation< dim, dim > &tria, const AdditionalData &additional_data=AdditionalData())
@ matrix
Contents is actually a matrix.
constexpr char T
constexpr types::blas_int zero
constexpr types::blas_int one
void mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
Definition l2.h:57
void partition(const SparsityPattern &sparsity_pattern, const unsigned int n_partitions, std::vector< unsigned int > &partition_indices, const Partitioner partitioner=Partitioner::metis)
void distribute_sparsity_pattern(DynamicSparsityPattern &dsp, const IndexSet &locally_owned_rows, const MPI_Comm mpi_comm, const IndexSet &locally_relevant_rows)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
Definition mpi.cc:105
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
Definition mpi.cc:120
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:466
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={})
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)
long double gamma(const unsigned int n)
int(& functions)(const void *v1, const void *v2)
void assemble(const MeshWorker::DoFInfoBox< dim, DOFINFO > &dinfo, A *assembler)
Definition loop.h:70
void refine_and_coarsen_fixed_fraction(::Triangulation< dim, spacedim > &tria, const ::Vector< Number > &criteria, const double top_fraction_of_error, const double bottom_fraction_of_error, const VectorTools::NormType norm_type=VectorTools::L1_norm)
::SolutionTransfer< dim, VectorType, spacedim > SolutionTransfer
void advance(std::tuple< I1, I2 > &t, const unsigned int n)
constexpr ProductType< Number, OtherNumber >::type scalar_product(const Tensor< rank, dim, Number > &left, const Tensor< rank, dim, OtherNumber > &right)
Definition tensor.h:2547