Reference documentation for deal.II version 9.6.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-4.h
Go to the documentation of this file.
1
393) const
394 *   {
395 *   double return_value = 0.0;
396 *   for (unsigned int i = 0; i < dim; ++i)
397 *   return_value += 4.0 * std::pow(p[i], 4.0);
398 *  
399 *   return return_value;
400 *   }
401 *  
402 *  
403 * @endcode
404 *
405 * As boundary values, we choose @f$x^2+y^2@f$ in 2d, and @f$x^2+y^2+z^2@f$ in 3d. This
406 * happens to be equal to the square of the vector from the origin to the
407 * point at which we would like to evaluate the function, irrespective of the
408 * dimension. So that is what we return:
409 *
410 * @code
411 *   template <int dim>
412 *   double BoundaryValues<dim>::value(const Point<dim> &p,
413 *   const unsigned int /*component*/) const
414 *   {
415 *   return p.square();
416 *   }
417 *  
418 *  
419 *  
420 * @endcode
421 *
422 *
423 * <a name="step_4-ImplementationofthecodeStep4codeclass"></a>
424 * <h3>Implementation of the <code>Step4</code> class</h3>
425 *
426
427 *
428 * Next for the implementation of the class template that makes use of the
429 * functions above. As before, we will write everything as templates that have
430 * a formal parameter <code>dim</code> that we assume unknown at the time we
431 * define the template functions. Only later, the compiler will find a
432 * declaration of <code>Step4@<2@></code> (in the <code>main</code> function,
433 * actually) and compile the entire class with <code>dim</code> replaced by 2,
434 * a process referred to as `instantiation of a template'. When doing so, it
435 * will also replace instances of <code>RightHandSide@<dim@></code> by
436 * <code>RightHandSide@<2@></code> and instantiate the latter class from the
437 * class template.
438 *
439
440 *
441 * In fact, the compiler will also find a declaration <code>Step4@<3@></code>
442 * in <code>main()</code>. This will cause it to again go back to the general
443 * <code>Step4@<dim@></code> template, replace all occurrences of
444 * <code>dim</code>, this time by 3, and compile the class a second time. Note
445 * that the two instantiations <code>Step4@<2@></code> and
446 * <code>Step4@<3@></code> are completely independent classes; their only
447 * common feature is that they are both instantiated from the same general
448 * template, but they are not convertible into each other, for example, and
449 * share no code (both instantiations are compiled completely independently).
450 *
451
452 *
453 *
454
455 *
456 *
457 * <a name="step_4-Step4Step4"></a>
458 * <h4>Step4::Step4</h4>
459 *
460
461 *
462 * After this introduction, here is the constructor of the <code>Step4</code>
463 * class. It specifies the desired polynomial degree of the finite elements
464 * and associates the DoFHandler to the triangulation just as in the previous
465 * example program, @ref step_3 "step-3":
466 *
467 * @code
468 *   template <int dim>
469 *   Step4<dim>::Step4()
470 *   : fe(/* polynomial degree = */ 1)
471 *   , dof_handler(triangulation)
472 *   {}
473 *  
474 *  
475 * @endcode
476 *
477 *
478 * <a name="step_4-Step4make_grid"></a>
479 * <h4>Step4::make_grid</h4>
480 *
481
482 *
483 * Grid creation is something inherently dimension dependent. However, as long
484 * as the domains are sufficiently similar in 2d or 3d, the library can
485 * abstract for you. In our case, we would like to again solve on the square
486 * @f$[-1,1]\times [-1,1]@f$ in 2d, or on the cube @f$[-1,1] \times [-1,1] \times
487 * [-1,1]@f$ in 3d; both can be termed GridGenerator::hyper_cube(), so we may
488 * use the same function in whatever dimension we are. Of course, the
489 * functions that create a hypercube in two and three dimensions are very much
490 * different, but that is something you need not care about. Let the library
491 * handle the difficult things.
492 *
493 * @code
494 *   template <int dim>
495 *   void Step4<dim>::make_grid()
496 *   {
497 *   GridGenerator::hyper_cube(triangulation, -1, 1);
498 *   triangulation.refine_global(4);
499 *  
500 *   std::cout << " Number of active cells: " << triangulation.n_active_cells()
501 *   << std::endl
502 *   << " Total number of cells: " << triangulation.n_cells()
503 *   << std::endl;
504 *   }
505 *  
506 * @endcode
507 *
508 *
509 * <a name="step_4-Step4setup_system"></a>
510 * <h4>Step4::setup_system</h4>
511 *
512
513 *
514 * This function looks exactly like in the previous example, although it
515 * performs actions that in their details are quite different if
516 * <code>dim</code> happens to be 3. The only significant difference from a
517 * user's perspective is the number of cells resulting, which is much higher
518 * in three than in two space dimensions!
519 *
520 * @code
521 *   template <int dim>
522 *   void Step4<dim>::setup_system()
523 *   {
524 *   dof_handler.distribute_dofs(fe);
525 *  
526 *   std::cout << " Number of degrees of freedom: " << dof_handler.n_dofs()
527 *   << std::endl;
528 *  
529 *   DynamicSparsityPattern dsp(dof_handler.n_dofs());
530 *   DoFTools::make_sparsity_pattern(dof_handler, dsp);
531 *   sparsity_pattern.copy_from(dsp);
532 *  
533 *   system_matrix.reinit(sparsity_pattern);
534 *  
535 *   solution.reinit(dof_handler.n_dofs());
536 *   system_rhs.reinit(dof_handler.n_dofs());
537 *   }
538 *  
539 *  
540 * @endcode
541 *
542 *
543 * <a name="step_4-Step4assemble_system"></a>
544 * <h4>Step4::assemble_system</h4>
545 *
546
547 *
548 * Unlike in the previous example, we would now like to use a non-constant
549 * right hand side function and non-zero boundary values. Both are tasks that
550 * are readily achieved with only a few new lines of code in the assemblage of
551 * the matrix and right hand side.
552 *
553
554 *
555 * More interesting, though, is the way we assemble matrix and right hand side
556 * vector dimension independently: there is simply no difference to the
557 * two-dimensional case. Since the important objects used in this function
558 * (quadrature formula, FEValues) depend on the dimension by way of a template
559 * parameter as well, they can take care of setting up properly everything for
560 * the dimension for which this function is compiled. By declaring all classes
561 * which might depend on the dimension using a template parameter, the library
562 * can make nearly all work for you and you don't have to care about most
563 * things.
564 *
565 * @code
566 *   template <int dim>
567 *   void Step4<dim>::assemble_system()
568 *   {
569 *   const QGauss<dim> quadrature_formula(fe.degree + 1);
570 *  
571 * @endcode
572 *
573 * We wanted to have a non-constant right hand side, so we use an object of
574 * the class declared above to generate the necessary data. Since this right
575 * hand side object is only used locally in the present function, we declare
576 * it here as a local variable:
577 *
578 * @code
579 *   RightHandSide<dim> right_hand_side;
580 *  
581 * @endcode
582 *
583 * Compared to the previous example, in order to evaluate the non-constant
584 * right hand side function we now also need the quadrature points on the
585 * cell we are presently on (previously, we only required values and
586 * gradients of the shape function from the FEValues object, as well as the
587 * quadrature weights, FEValues::JxW() ). We can tell the FEValues object to
588 * do for us by also giving it the #update_quadrature_points flag:
589 *
590 * @code
591 *   FEValues<dim> fe_values(fe,
592 *   quadrature_formula,
593 *   update_values | update_gradients |
594 *   update_quadrature_points | update_JxW_values);
595 *  
596 * @endcode
597 *
598 * We then again define the same abbreviation as in the previous program.
599 * The value of this variable of course depends on the dimension which we
600 * are presently using, but the FiniteElement class does all the necessary
601 * work for you and you don't have to care about the dimension dependent
602 * parts:
603 *
604 * @code
605 *   const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
606 *  
607 *   FullMatrix<double> cell_matrix(dofs_per_cell, dofs_per_cell);
608 *   Vector<double> cell_rhs(dofs_per_cell);
609 *  
610 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
611 *  
612 * @endcode
613 *
614 * Next, we again have to loop over all cells and assemble local
615 * contributions. Note, that a cell is a quadrilateral in two space
616 * dimensions, but a hexahedron in 3d. In fact, the
617 * <code>active_cell_iterator</code> data type is something different,
618 * depending on the dimension we are in, but to the outside world they look
619 * alike and you will probably never see a difference. In any case, the real
620 * type is hidden by using `auto`:
621 *
622 * @code
623 *   for (const auto &cell : dof_handler.active_cell_iterators())
624 *   {
625 *   fe_values.reinit(cell);
626 *  
627 *   cell_matrix = 0;
628 *   cell_rhs = 0;
629 *  
630 * @endcode
631 *
632 * Now we have to assemble the local matrix and right hand side. This is
633 * done exactly like in the previous example, but now we revert the
634 * order of the loops (which we can safely do since they are independent
635 * of each other) and merge the loops for the local matrix and the local
636 * vector as far as possible to make things a bit faster.
637 *
638
639 *
640 * Assembling the right hand side presents the only significant
641 * difference to how we did things in @ref step_3 "step-3": Instead of using a
642 * constant right hand side with value 1, we use the object representing
643 * the right hand side and evaluate it at the quadrature points:
644 *
645 * @code
646 *   for (const unsigned int q_index : fe_values.quadrature_point_indices())
647 *   for (const unsigned int i : fe_values.dof_indices())
648 *   {
649 *   for (const unsigned int j : fe_values.dof_indices())
650 *   cell_matrix(i, j) +=
651 *   (fe_values.shape_grad(i, q_index) * // grad phi_i(x_q)
652 *   fe_values.shape_grad(j, q_index) * // grad phi_j(x_q)
653 *   fe_values.JxW(q_index)); // dx
654 *  
655 *   const auto &x_q = fe_values.quadrature_point(q_index);
656 *   cell_rhs(i) += (fe_values.shape_value(i, q_index) * // phi_i(x_q)
657 *   right_hand_side.value(x_q) * // f(x_q)
658 *   fe_values.JxW(q_index)); // dx
659 *   }
660 * @endcode
661 *
662 * As a final remark to these loops: when we assemble the local
663 * contributions into <code>cell_matrix(i,j)</code>, we have to multiply
664 * the gradients of shape functions @f$i@f$ and @f$j@f$ at point number
665 * q_index and
666 * multiply it with the scalar weights JxW. This is what actually
667 * happens: <code>fe_values.shape_grad(i,q_index)</code> returns a
668 * <code>dim</code> dimensional vector, represented by a
669 * <code>Tensor@<1,dim@></code> object, and the operator* that
670 * multiplies it with the result of
671 * <code>fe_values.shape_grad(j,q_index)</code> makes sure that the
672 * <code>dim</code> components of the two vectors are properly
673 * contracted, and the result is a scalar floating point number that
674 * then is multiplied with the weights. Internally, this operator* makes
675 * sure that this happens correctly for all <code>dim</code> components
676 * of the vectors, whether <code>dim</code> be 2, 3, or any other space
677 * dimension; from a user's perspective, this is not something worth
678 * bothering with, however, making things a lot simpler if one wants to
679 * write code dimension independently.
680 *
681
682 *
683 * With the local systems assembled, the transfer into the global matrix
684 * and right hand side is done exactly as before, but here we have again
685 * merged some loops for efficiency:
686 *
687 * @code
688 *   cell->get_dof_indices(local_dof_indices);
689 *   for (const unsigned int i : fe_values.dof_indices())
690 *   {
691 *   for (const unsigned int j : fe_values.dof_indices())
692 *   system_matrix.add(local_dof_indices[i],
693 *   local_dof_indices[j],
694 *   cell_matrix(i, j));
695 *  
696 *   system_rhs(local_dof_indices[i]) += cell_rhs(i);
697 *   }
698 *   }
699 *  
700 * @endcode
701 *
702 * As the final step in this function, we wanted to have non-homogeneous
703 * boundary values in this example, unlike the one before. This is a simple
704 * task, we only have to replace the Functions::ZeroFunction used there by an
705 * object of the class which describes the boundary values we would like to
706 * use (i.e. the <code>BoundaryValues</code> class declared above):
707 *
708
709 *
710 * The function VectorTools::interpolate_boundary_values() will only work
711 * on faces that have been marked with boundary indicator 0 (because that's
712 * what we say the function should work on with the second argument below).
713 * If there are faces with boundary id other than 0, then the function
714 * interpolate_boundary_values() will do nothing on these faces. For
715 * the Laplace equation doing nothing is equivalent to assuming that
716 * on those parts of the boundary a zero Neumann boundary condition holds.
717 *
718 * @code
719 *   std::map<types::global_dof_index, double> boundary_values;
721 *   types::boundary_id(0),
722 *   BoundaryValues<dim>(),
723 *   boundary_values);
724 *   MatrixTools::apply_boundary_values(boundary_values,
725 *   system_matrix,
726 *   solution,
727 *   system_rhs);
728 *   }
729 *  
730 *  
731 * @endcode
732 *
733 *
734 * <a name="step_4-Step4solve"></a>
735 * <h4>Step4::solve</h4>
736 *
737
738 *
739 * Solving the linear system of equations is something that looks almost
740 * identical in most programs. In particular, it is dimension independent, so
741 * this function is copied verbatim from the previous example.
742 *
743 * @code
744 *   template <int dim>
745 *   void Step4<dim>::solve()
746 *   {
747 *   SolverControl solver_control(1000, 1e-6 * system_rhs.l2_norm());
748 *   SolverCG<Vector<double>> solver(solver_control);
749 *   solver.solve(system_matrix, solution, system_rhs, PreconditionIdentity());
750 *  
751 *   std::cout << " " << solver_control.last_step()
752 *   << " CG iterations needed to obtain convergence." << std::endl;
753 *   }
754 *  
755 *  
756 * @endcode
757 *
758 *
759 * <a name="step_4-Step4output_results"></a>
760 * <h4>Step4::output_results</h4>
761 *
762
763 *
764 * This function also does what the respective one did in @ref step_3 "step-3". No changes
765 * here for dimension independence either.
766 *
767
768 *
769 * Since the program will run both 2d and 3d versions of the Laplace solver,
770 * we use the dimension in the filename to generate distinct filenames for
771 * each run (in a better program, one would check whether <code>dim</code> can
772 * have other values than 2 or 3, but we neglect this here for the sake of
773 * brevity).
774 *
775 * @code
776 *   template <int dim>
777 *   void Step4<dim>::output_results() const
778 *   {
779 *   DataOut<dim> data_out;
780 *  
781 *   data_out.attach_dof_handler(dof_handler);
782 *   data_out.add_data_vector(solution, "solution");
783 *  
784 *   data_out.build_patches();
785 *  
786 *   std::ofstream output(dim == 2 ? "solution-2d.vtk" : "solution-3d.vtk");
787 *   data_out.write_vtk(output);
788 *   }
789 *  
790 *  
791 *  
792 * @endcode
793 *
794 *
795 * <a name="step_4-Step4run"></a>
796 * <h4>Step4::run</h4>
797 *
798
799 *
800 * This is the function which has the top-level control over everything. Apart
801 * from one line of additional output, it is the same as for the previous
802 * example.
803 *
804 * @code
805 *   template <int dim>
806 *   void Step4<dim>::run()
807 *   {
808 *   std::cout << "Solving problem in " << dim << " space dimensions."
809 *   << std::endl;
810 *  
811 *   make_grid();
812 *   setup_system();
813 *   assemble_system();
814 *   solve();
815 *   output_results();
816 *   }
817 *  
818 *  
819 * @endcode
820 *
821 *
822 * <a name="step_4-Thecodemaincodefunction"></a>
823 * <h3>The <code>main</code> function</h3>
824 *
825
826 *
827 * And this is the main function. It also looks mostly like in @ref step_3 "step-3", but if
828 * you look at the code below, note how we first create a variable of type
829 * <code>Step4@<2@></code> (forcing the compiler to compile the class template
830 * with <code>dim</code> replaced by <code>2</code>) and run a 2d simulation,
831 * and then we do the whole thing over in 3d.
832 *
833
834 *
835 * In practice, this is probably not what you would do very frequently (you
836 * probably either want to solve a 2d problem, or one in 3d, but not both at
837 * the same time). However, it demonstrates the mechanism by which we can
838 * simply change which dimension we want in a single place, and thereby force
839 * the compiler to recompile the dimension independent class templates for the
840 * dimension we request. The emphasis here lies on the fact that we only need
841 * to change a single place. This makes it rather trivial to debug the program
842 * in 2d where computations are fast, and then switch a single place to a 3 to
843 * run the much more computing intensive program in 3d for `real'
844 * computations.
845 *
846
847 *
848 * Each of the two blocks is enclosed in braces to make sure that the
849 * <code>laplace_problem_2d</code> variable goes out of scope (and releases
850 * the memory it holds) before we move on to allocate memory for the 3d
851 * case. Without the additional braces, the <code>laplace_problem_2d</code>
852 * variable would only be destroyed at the end of the function, i.e. after
853 * running the 3d problem, and would needlessly hog memory while the 3d run
854 * could actually use it.
855 *
856 * @code
857 *   int main()
858 *   {
859 *   {
860 *   Step4<2> laplace_problem_2d;
861 *   laplace_problem_2d.run();
862 *   }
863 *  
864 *   {
865 *   Step4<3> laplace_problem_3d;
866 *   laplace_problem_3d.run();
867 *   }
868 *  
869 *   return 0;
870 *   }
871 * @endcode
872<a name="step_4-Results"></a><h1>Results</h1>
873
874
875
876The output of the program looks as follows (the number of iterations
877may vary by one or two, depending on your computer, since this is
878often dependent on the round-off accuracy of floating point
879operations, which differs between processors):
880@code
881Solving problem in 2 space dimensions.
882 Number of active cells: 256
883 Total number of cells: 341
884 Number of degrees of freedom: 289
885 19 CG iterations needed to obtain convergence.
886Solving problem in 3 space dimensions.
887 Number of active cells: 4096
888 Total number of cells: 4681
889 Number of degrees of freedom: 4913
890 20 CG iterations needed to obtain convergence.
891@endcode
892It is obvious that in three spatial dimensions the number of cells and
893therefore also the number of degrees of freedom is
894much higher. What cannot be seen here, is that besides this higher
895number of rows and columns in the matrix, there are also significantly
896more entries per row of the matrix in three space
897dimensions. Together, this leads to a much higher numerical effort for
898solving the system of equation, which you can feel in the run time of the two
899solution steps when you actually run the program.
900
901
902
903The program produces two files: <code>solution-2d.vtk</code> and
904<code>solution-3d.vtk</code>, which can be viewed using the programs
905VisIt or Paraview (in case you do not have these programs, you can easily
906change the
907output format in the program to something which you can view more
908easily). Visualizing solutions is a bit of an art, but it can also be fun, so
909you should play around with your favorite visualization tool to get familiar
910with its functionality. Here's what I have come up with for the 2d solution:
911
912<p align="center">
913 <img src="https://www.dealii.org/images/steps/developer/step-4.solution-2d.png" alt="">
914</p>
915
916(See also <a href="https://www.math.colostate.edu/~bangerth/videos.676.11.html">video lecture 11</a>, <a href="https://www.math.colostate.edu/~bangerth/videos.676.32.html">video lecture 32</a>.)
917The picture shows the solution of the problem under consideration as
918a 3D plot. As can be seen, the solution is almost flat in the interior
919of the domain and has a higher curvature near the boundary. This, of
920course, is due to the fact that for Laplace's equation the curvature
921of the solution is equal to the right hand side and that was chosen as
922a quartic polynomial which is nearly zero in the interior and is only
923rising sharply when approaching the boundaries of the domain; the
924maximal values of the right hand side function are at the corners of
925the domain, where also the solution is moving most rapidly.
926It is also nice to see that the solution follows the desired quadratic
927boundary values along the boundaries of the domain.
928It can also be useful to verify a computed solution against an analytical
929solution. For an explanation of this technique, see @ref step_7 "step-7".
930
931On the other hand, even though the picture does not show the mesh lines
932explicitly, you can see them as little kinks in the solution. This clearly
933indicates that the solution hasn't been computed to very high accuracy and
934that to get a better solution, we may have to compute on a finer mesh.
935
936In three spatial dimensions, visualization is a bit more difficult. The left
937picture shows the solution and the mesh it was computed on on the surface of
938the domain. This is nice, but it has the drawback that it completely hides
939what is happening on the inside. The picture on the right is an attempt at
940visualizing the interior as well, by showing surfaces where the solution has
941constant values (as indicated by the legend at the top left). Isosurface
942pictures look best if one makes the individual surfaces slightly transparent
943so that it is possible to see through them and see what's behind.
944
945<table width="60%" align="center">
946 <tr>
947 <td align="center">
948 <img src="https://www.dealii.org/images/steps/developer/step-4.solution-3d.png" alt="">
949 </td>
950 <td align="center">
951 <img src="https://www.dealii.org/images/steps/developer/step-4.contours-3d.png" alt="">
952 </td>
953 </tr>
954</table>
955
956@note
957A final remark on visualization: the idea of visualization is to give insight,
958which is not the same as displaying information. In particular, it is easy to
959overload a picture with information, but while it shows more information it
960makes it also more difficult to glean insight. As an example, the program I
961used to generate these pictures, VisIt, by default puts tick marks on every
962axis, puts a big fat label "X Axis" on the @f$x@f$ axis and similar for the other
963axes, shows the file name from which the data was taken in the top left and
964the name of the user doing so and the time and date on the bottom right. None
965of this is important
966here: the axes are equally easy to make out because the tripod at the bottom
967left is still visible, and we know from the program that the domain is
968@f$[-1,1]^3@f$, so there is no need for tick marks. As a consequence, I have
969switched off all the extraneous stuff in the picture: the art of visualization
970is to reduce the picture to those parts that are important to see what one
971wants to see, but no more.
972
973
974
975<a name="step-4-postprocessing"></a>
976<a name="step_4-PostprocessingWhattodowiththesolution"></a><h3> Postprocessing: What to do with the solution? </h3>
977
978
979This tutorial -- like most of the other programs -- principally only shows how
980to numerically approximate the solution of a partial differential equation, and
981then how to visualize this solution graphically. But
982solving a PDE is of course not the goal in most practical applications (unless
983you are a numerical methods developer and the *method* is the goal): We generally
984want to solve a PDE because we want to *extract information* from it. Examples
985for what people are interested in from solutions include the following:
986- Let's say you solve the equations of elasticity (which we will do in @ref step_8 "step-8"),
987 then that's presumably because you want to know about the deformation of an
988 elastic object under a given load. From an engineering perspective, what you
989 then presumably want to learn is the degree of deformation of the object,
990 say at a specific point; or you may want to know the maximum
991 [stress](https://en.wikipedia.org/wiki/Stress_(mechanics)) in order to
992 determine whether the applied load exceeds the safe maximal stress the
993 material can withstand.
994- If you are solving fluid flow problems (such as in @ref step_22 "step-22", @ref step_57 "step-57", @ref step_67 "step-67",
995 and @ref step_69 "step-69"), then you might be interested in the fluid velocity at specific
996 points, and oftentimes the forces the fluid exerts on the boundary of the
997 fluid domain. The latter is important in many applications: If the fluid
998 in question is the air flowing around an airplane, then we are typically
999 interested in the drag and lift forces on the fuselage and wings. If the
1000 fluid is water flowing around a ship, then we typically care about
1001 the drag force on the ship.
1002- If you are solving the Maxwell equations of electromagnetics, you are
1003 typically interested in how much energy is radiated in certain directions
1004 (say, in order to know the range of information transmission via an
1005 antenna, or to determine the
1006 [radar cross section](https://en.wikipedia.org/wiki/Radar_cross_section)
1007 of planes or ships).
1008
1009The point here is that from an engineering perspective, *solving* the
1010PDE is only the first step. The second step is to *evaluate* the computed
1011solution in order to extract relevant numbers that allow us to
1012either *optimize a design*, or to *make decisions*. This second step is often
1013called "postprocessing the solution".
1014
1015This program does not solve a solid or fluid mechanics problem, so we
1016should try to illustrate postprocessing with something that makes sense
1017in the context of the equation we solve here. The Poisson equation in two
1018space dimensions is a model for the vertical deformation of a membrane
1019that is clamped at the boundary and is subject to a vertical force. For this
1020kind of situation, it makes sense to evaluate the *average vertical
1021displacement*,
1022@f[
1023 \bar u_h = \frac{\int_\Omega u_h(\mathbf x) \, dx}{|\Omega|},
1024@f]
1025where @f$|\Omega| = \int_\Omega 1 \, dx@f$ is the area of the domain. To compute
1026@f$\bar u_h@f$, as usual we replace integrals over the domain by a sum of integrals
1027over cells,
1028@f[
1029 \int_\Omega u_h(\mathbf x) \, dx
1030 =
1031 \sum_K \int_K u_h(\mathbf x) \, dx,
1032@f]
1033and then integrals over cells are approximated by quadrature:
1034@f{align*}{
1035 \int_\Omega u_h(\mathbf x) \, dx
1036 &=
1037 \sum_K \int_K u_h(\mathbf x) \, dx,
1038 \\
1039 &=
1040 \sum_K \sum_q u_h(\mathbf x_q^K) w_q^K,
1041@f}
1042where @f$w_q^K@f$ is the weight of the @f$q@f$th quadrature point evaluated on
1043cell @f$K@f$. All of this is as always provided by the FEValues class -- the
1044entry point for all integrals in deal.II.
1045
1046The actual implementation of this is straightforward once you know how
1047to get the values of the solution @f$u@f$ at the quadrature points of a cell.
1048This functionality is provided by FEValues::get_function_values(), a
1049function that takes a global vector of nodal values as input and returns
1050a vector of function values at the quadrature points of the current cell.
1051Using this function, to see how it all works together you can
1052place the following code snippet anywhere in the program after the
1053solution has been computed (the `output_results()` function seems like a good
1054place to also do postprocessing, for example):
1055@code
1056 QGauss<dim> quadrature_formula(fe.degree + 1);
1057 FEValues<dim> fe_values(fe,
1058 quadrature_formula,
1060
1061 std::vector<double> solution_values(quadrature_formula.size());
1062 double integral_of_u = 0;
1063 double volume_of_omega = 0;
1064
1065 for (const auto &cell : dof_handler.active_cell_iterators())
1066 {
1067 fe_values.reinit(cell);
1068 fe_values.get_function_values(solution, solution_values);
1069
1070 for (const unsigned int q_point : fe_values.quadrature_point_indices())
1071 {
1072 integral_of_u += solution_values[q_point] * fe_values.JxW(q_point);
1073 volume_of_omega += 1 * fe_values.JxW(q_point);
1074 }
1075 }
1076 std::cout << " Mean value of u=" << integral_of_u / volume_of_omega
1077 << std::endl;
1078@endcode
1079In this code snippet, we also compute the volume (or, since we are currently
1080thinking about a two-dimensional situation: the area) @f$|\Omega|@f$ by computing
1081the integral @f$|\Omega| = \int_\Omega 1 \, dx@f$ in exactly the same way, via
1082quadrature. (We could avoid having to compute @f$|\Omega|@f$ by hand here, using the
1083fact that deal.II has a function for this: GridTools::volume(). That said,
1084it is efficient to compute the two integrals
1085concurrently in the same loop, and so that's what we do.)
1086
1087This program of course also solves the same Poisson equation in three space
1088dimensions. In this situation, the Poisson equation is often used as a model
1089for diffusion of either a physical species (say, of ink in a tank of water,
1090or a pollutant in the air) or of energy (specifically, of thermal energy in
1091a solid body). In that context, the quantity
1092@f[
1093 \Phi_h = \int_{\partial\Omega} \nabla u_h(\mathbf x) \cdot \mathbf n(\mathbf x) \; dx
1094@f]
1095is the *flux* of this species or energy across the boundary. (In actual
1096physical models, one would also have to multiply the right hand side by
1097a diffusivity or conductivity constant, but let us ignore this here.) In
1098much the same way as before, we compute such integrals by splitting
1099it over integrals of *faces* of cells, and then applying quadrature:
1100@f{align*}{
1101 \Phi_h
1102 &=
1103 \int_{\partial\Omega} \nabla u_h(\mathbf x) \cdot \mathbf n(\mathbf x) \; dx
1104 \\
1105 &=
1106 \sum_K
1107 \sum_{f \in \text{faces of @f$K@f$}, f\subset\partial\Omega}
1108 \int_f \nabla u_h(\mathbf x) \cdot \mathbf n(\mathbf x) \; dx
1109 \\
1110 &=
1111 \sum_K
1112 \sum_{f \in \text{faces of @f$K@f$}, f\subset\partial\Omega}
1113 \sum_q \nabla u_h(\mathbf x_q^f) \cdot \mathbf n(\mathbf x_q^f) w_q^f,
1114@f}
1115where now @f$\mathbf x_q^f@f$ are the quadrature points located on face @f$f@f$,
1116and @f$w_q^f@f$ are the weights associated with these faces. The second
1117of the sum symbols loops over all faces of cell @f$K@f$, but restricted to
1118those that are actually at the boundary.
1119
1120This all is easily implemented by the following code that replaces the use of the
1121FEValues class (which is used for integrating over cells -- i.e., domain integrals)
1122by the FEFaceValues class (which is used for integrating over faces -- i.e.,
1123boundary integrals):
1124@code
1125 QGauss<dim - 1> face_quadrature_formula(fe.degree + 1);
1126 FEFaceValues<dim> fe_face_values(fe,
1127 face_quadrature_formula,
1128 update_gradients | update_normal_vectors |
1129 update_JxW_values);
1130
1131 std::vector<Tensor<1, dim>> solution_gradients(face_quadrature_formula.size());
1132 double flux = 0;
1133
1134 for (const auto &cell : dof_handler.active_cell_iterators())
1135 for (const auto &face : cell->face_iterators())
1136 if (face->at_boundary())
1137 {
1138 fe_face_values.reinit(cell, face);
1139 fe_face_values.get_function_gradients(solution, solution_gradients);
1140
1141 for (const unsigned int q_point :
1142 fe_face_values.quadrature_point_indices())
1143 {
1144 flux += solution_gradients[q_point] *
1145 fe_face_values.normal_vector(q_point) *
1146 fe_face_values.JxW(q_point);
1147 }
1148 }
1149 std::cout << " Flux=" << flux << std::endl;
1150@endcode
1151
1152If you add these two code snippets to the code, you will get output like the
1153following when you run the program:
1154@code
1155Solving problem in 2 space dimensions.
1156 Number of active cells: 256
1157 Total number of cells: 341
1158 Number of degrees of freedom: 289
1159 26 CG iterations needed to obtain convergence.
1160 Mean value of u=1.33303
1161 Flux=-3.68956
1162Solving problem in 3 space dimensions.
1163 Number of active cells: 4096
1164 Total number of cells: 4681
1165 Number of degrees of freedom: 4913
1166 30 CG iterations needed to obtain convergence.
1167 Mean value of u=1.58058
1168 Flux=-8.29435
1169@endcode
1170
1171This makes some sense: If you look, for example, at the 2d output above,
1172the solution varies between values of 1 and 2, but with a larger part of the
1173solution closer to one than two; so an average value of 1.33 for the mean value
1174is reasonable. For the flux, recall that @f$\nabla u \cdot \mathbf n@f$ is the
1175directional derivative in the normal direction -- in other words, how the
1176solution changes as we move from the interior of the domain towards the
1177boundary. If you look at the 2d solution, you will realize that for most parts
1178of the boundary, the solution *decreases* as we approach the boundary, so the
1179normal derivative is negative -- so if we integrate along the boundary, we
1180should expect (and obtain!) a negative value.
1181
1182
1183
1184<a name="step-4-extensions"></a>
1185<a name="step_4-Possibilitiesforextensions"></a><h3>Possibilities for extensions</h3>
1186
1187
1188There are many ways with which one can play with this program. The simpler
1189ones include essentially all the possibilities already discussed in the
1190<a href="step_3.html#extensions" target="body">Possibilities for extensions in the documentation of step 3</a>,
1191except that you will have to think about whether something now also applies
1192to the 3d case discussed in the current program.
1193
1194It is also worthwhile considering the postprocessing options discussed
1195above. The documentation states two numbers (the mean value and the
1196normal flux) for both the 2d and 3d cases. Can we trust these
1197numbers? We have convinced ourselves that at least the mean value
1198is reasonable, and that the sign of the flux is probably correct.
1199But are these numbers accurate?
1200
1201A general rule is that we should never trust a number unless we have
1202verified it in some way. From the theory of finite element methods,
1203we know that as we make the mesh finer and finer, the numerical
1204solution @f$u_h@f$ we compute here must converge to the exact solution
1205@f$u@f$. As a consequence, we also expect that @f$\bar u_h \rightarrow \bar u@f$
1206and @f$\Phi_h \rightarrow \Phi@f$, but that does not mean that for any
1207given mesh @f$\bar u_h@f$ or @f$\Phi_h@f$ are particularly accurate approximations.
1208
1209To test this kind of thing, we have already considered the convergence of
1210a point value in @ref step_3 "step-3". We can do the same here by selecting how many
1211times the mesh is globally refined in the `make_grid()` function of this
1212program. For the mean value of the solution, we then get the following
1213numbers:
1214 <table align="center" class="doxtable">
1215 <tr> <th># of refinements</th>
1216 <th>@f$\bar u_h@f$ in 2d</th>
1217 <th>@f$\bar u_h@f$ in 3d</th>
1218 </tr>
1219 <tr> <td>4</td> <td>1.33303</td> <td>1.58058</td> </tr>
1220 <tr> <td>5</td> <td>1.33276</td> <td>1.57947</td> </tr>
1221 <tr> <td>6</td> <td>1.3327</td> <td>1.5792</td> </tr>
1222 <tr> <td>7</td> <td>1.33269</td> <td>1.57914</td> </tr>
1223 <tr> <td>8</td> <td>1.33268</td> <td></td> </tr>
1224 <tr> <td>9</td> <td>1.33268</td> <td></td> </tr>
1225 </table>
1226I did not have the patience to run the last two values for the 3d case --
1227one needs quite a fine mesh for this, with correspondingly long run times.
1228But we can be reasonably assured that values around 1.33 (for the 2d case)
1229and 1.58 (for the 3d case) are about right -- and at least for engineering
1230applications, three digits of accuracy are good enough.
1231
1232The situation looks very different for the flux. Here, we get results
1233such as the following:
1234 <table align="center" class="doxtable">
1235 <tr> <th># of refinements</th>
1236 <th>@f$\Phi_h@f$ in 2d</th>
1237 <th>@f$\Phi_h@f$ in 3d</th>
1238 </tr>
1239 <tr> <td>4</td> <td>-3.68956</td> <td>-8.29435</td> </tr>
1240 <tr> <td>5</td> <td>-4.90147</td> <td>-13.0691</td> </tr>
1241 <tr> <td>6</td> <td>-5.60745</td> <td>-15.9171</td> </tr>
1242 <tr> <td>7</td> <td>-5.99111</td> <td>-17.4918</td> </tr>
1243 <tr> <td>8</td> <td>-6.19196</td> <td></td> </tr>
1244 <tr> <td>9</td> <td>-6.29497</td> <td></td> </tr>
1245 <tr> <td>10</td> <td>-6.34721</td> <td></td> </tr>
1246 <tr> <td>11</td> <td>-6.37353</td> <td></td> </tr>
1247 </table>
1248So this is not great. For the 2d case, we might infer that perhaps
1249a value around -6.4 might be right if we just refine the mesh enough --
1250though 11 refinements already leads to some 4,194,304 cells. In any
1251case, the first number (the one shown in the beginning where we
1252discussed postprocessing) was off by almost a factor of 2!
1253
1254For the 3d case, the last number shown was on a mesh with 2,097,152
1255cells; the next one would have had 8 times as many cells. In any case, the
1256numbers mean that we can't even be sure
1257that the first digit of that last number is correct! In other words,
1258it was worth checking, or we would have just believed all of these
1259numbers. In fact, that last column isn't even doing a particularly
1260good job convincing us that the code might be correctly implemented.
1261
1262If you keep reading through the other tutorial programs, you will find many ways
1263to make these sorts of computations more accurate and to come to
1264believe that the flux actually does converge to its correct value.
1265For example, we can dramatically increase the accuracy of the computation
1266by using adaptive mesh refinement (@ref step_6 "step-6") near the boundary, and
1267in particular by using higher polynomial degree finite elements (also
1268@ref step_6 "step-6", but also @ref step_7 "step-7"). Using the latter, using cubic elements
1269(polynomial degree 3), we can actually compute the flux pretty
1270accurately even in 3d: @f$\Phi_h=-19.0148@f$ with 4 global refinement steps,
1271and @f$\Phi_h=-19.1533@f$ with 5 refinement steps. These numbers are already
1272pretty close together and give us a reasonable idea of the first
1273two correct digits of the "true" answer.
1274
1275@note We would be remiss to not also comment on the fact that there
1276 are good theoretical reasons why computing the flux accurately
1277 appears to be so much more difficult than the average value.
1278 This has to do with the fact that finite element theory
1279 provides us with the estimate
1280 @f$\|u-u_h\|_{L_2(\Omega)} \le C h^2 \|\nabla^2u\|_{L_2(\Omega)}@f$
1281 when using the linear elements this program uses -- that is, for
1282 every global mesh refinement, @f$h@f$ is reduced by a factor of two
1283 and the error goes down by a factor of 4. Now, the @f$L_2@f$ error is
1284 not equivalent to the error in the mean value, but the two are
1285 related: They are both integrals over the domain, using the *value*
1286 of the solution. We expect the mean value to converge no worse than
1287 the @f$L_2@f$ norm of the error. At the same time, theory also provides
1288 us with this estimate:
1289 @f$\|\nabla (u-u_h)\|_{L_2(\partial\Omega)} \le
1290 C h^{1/2} \|\nabla^2u\|_{L_2(\Omega)}@f$. The move from values to
1291 gradients reduces the convergence rates by one order, and the move
1292 from domain to boundary by another half order. Here, then, each
1293 refinement step reduces the error not by a factor of 4 any more,
1294 by only by a factor of @f$\sqrt{2} \approx 1.4@f$. It takes a lot
1295 of global refinement steps to reduce the error by, say, a factor
1296 ten or hundred, and this is reflected in the very slow convergence
1297 evidenced by the table. On the other hand, for cubic elements (i.e.,
1298 polynomial degree 3), we would get
1299 @f$\|u-u_h\|_{L_2(\Omega)} \le C h^4 \|\nabla^4u\|_{L_2(\Omega)}@f$
1300 and after reduction by 1.5 orders, we would still have
1301 @f$\|\nabla (u-u_h)\|_{L_2(\partial\Omega)} \le
1302 C h^{2+1/2} \|\nabla^4u\|_{L_2(\Omega)}@f$. This rate,
1303 @f${\cal O}(h^{2.5})@f$ is still quite rapid, and it is perhaps not
1304 surprising that we get much better answers with these higher
1305 order elements. This also illustrates that when trying to
1306 approximate anything that relates to a gradient of the solution,
1307 using linear elements (polynomial degree one) is really not a
1308 good choice at all.
1309
1310@note In this very specific case, it turns out that we can actually
1311 compute the exact value of @f$\Phi@f$. This is because for the Poisson
1312 equation we compute the solution of here, @f$-\Delta u = f@f$, we can
1313 integrate over the domain, @f$-\int_\Omega \Delta u = \int_\Omega f@f$,
1314 and then use that @f$\Delta = \text{div}\;\text{grad}@f$; this allows
1315 us to use the
1316 divergence theorem followed by multiplying by minus one to find
1317 @f$\int_{\partial\Omega} \nabla u \cdot n = -\int_\Omega f@f$. The
1318 left hand side happens to be @f$\Phi@f$. For the specific right
1319 hand side @f$f(x_1,x_2)=4(x_1^4+x_2^4)@f$ we use in 2d, we then
1320 get @f$-\int_\Omega f = -\int_{-1}^{1} \int_{-1}^{1} 4(x_1^4+x_2^4) \; dx_2\; dx_1
1321 = -16 \left[\int_{-1}^{1} x^4 \; dx\right] = -16\times\frac 25@f$,
1322 which has a numerical value of exactly -6.4 -- right on with our
1323 guess above. In 3d, we can do the same and get that the exact
1324 value is
1325 @f$-\int_\Omega f =
1326 -\int_{-1}^{1} \int_{-1}^{1} \int_{-1}^{1} 4(x_1^4+x_2^4+x_3^4) \; dx_3 \; dx_2\; dx_1
1327 = -48\times\frac 25=-19.2@f$. What we found with cubic elements
1328 is then quite close to this exact value. Of course, in practice
1329 we almost never have exact values to compare with: If we could
1330 compute something on a piece of paper, we wouldn't have to solve
1331 the PDE numerically. But these sorts of situations make for excellent
1332 test cases that help us verify that our numerical solver works
1333 correctly. In many other cases, the literature contains
1334 numbers where others have already computed an answer accurately
1335 using their own software, and these are also often useful to
1336 compare against in verifying the correctness of our codes.
1337 *
1338 *
1339<a name="step_4-PlainProg"></a>
1340<h1> The plain program</h1>
1341@include "step-4.cc"
1342*/
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
Definition point.h:111
constexpr numbers::NumberTraits< Number >::real_type square() const
Point< 3 > center
Point< 2 > second
Definition grid_out.cc:4624
Point< 2 > first
Definition grid_out.cc:4623
unsigned int level
Definition grid_out.cc:4626
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
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_JxW_values
Transformed quadrature weights.
CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt K
double volume(const Triangulation< dim, spacedim > &tria)
@ matrix
Contents is actually a matrix.
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:74
void apply_boundary_values(const std::map< types::global_dof_index, number > &boundary_values, SparseMatrix< number > &matrix, Vector< number > &solution, Vector< number > &right_hand_side, const bool eliminate_columns=true)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:191
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
T sum(const T &t, const MPI_Comm mpi_communicator)
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)
DEAL_II_HOST constexpr TableIndices< 2 > merge(const TableIndices< 2 > &previous_indices, const unsigned int new_index, const unsigned int position)
int(& functions)(const void *v1, const void *v2)
void assemble(const MeshWorker::DoFInfoBox< dim, DOFINFO > &dinfo, A *assembler)
Definition loop.h:70
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)