deal.II version GIT relicensing-2173-gae8fc9d14b 2024-11-24 06:40:00+00:00
\(\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-5.h
Go to the documentation of this file.
1 1)
402 *   , dof_handler(triangulation)
403 *   {}
404 *  
405 *  
406 *  
407 * @endcode
408 *
409 *
410 * <a name="step_5-Step5setup_system"></a>
411 * <h4>Step5::setup_system</h4>
412 *
413
414 *
415 * This is the function <code>make_grid</code> from the previous
416 * example, minus the generation of the grid. Everything else is unchanged:
417 *
418 * @code
419 *   template <int dim>
420 *   void Step5<dim>::setup_system()
421 *   {
422 *   dof_handler.distribute_dofs(fe);
423 *  
424 *   std::cout << " Number of degrees of freedom: " << dof_handler.n_dofs()
425 *   << std::endl;
426 *  
427 *   DynamicSparsityPattern dsp(dof_handler.n_dofs());
428 *   DoFTools::make_sparsity_pattern(dof_handler, dsp);
429 *   sparsity_pattern.copy_from(dsp);
430 *  
431 *   system_matrix.reinit(sparsity_pattern);
432 *  
433 *   solution.reinit(dof_handler.n_dofs());
434 *   system_rhs.reinit(dof_handler.n_dofs());
435 *   }
436 *  
437 *  
438 *  
439 * @endcode
440 *
441 *
442 * <a name="step_5-Step5assemble_system"></a>
443 * <h4>Step5::assemble_system</h4>
444 *
445
446 *
447 * As in the previous examples, this function is not changed much with regard
448 * to its functionality, but there are still some optimizations which we will
449 * show. For this, it is important to note that if efficient solvers are used
450 * (such as the preconditioned CG method), assembling the matrix and right hand
451 * side can take a comparable time, and you should think about using one or
452 * two optimizations at some places.
453 *
454
455 *
456 * The first parts of the function are completely unchanged from before:
457 *
458 * @code
459 *   template <int dim>
460 *   void Step5<dim>::assemble_system()
461 *   {
462 *   const QGauss<dim> quadrature_formula(fe.degree + 1);
463 *  
464 *   FEValues<dim> fe_values(fe,
465 *   quadrature_formula,
468 *  
469 *   const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
470 *  
471 *   FullMatrix<double> cell_matrix(dofs_per_cell, dofs_per_cell);
472 *   Vector<double> cell_rhs(dofs_per_cell);
473 *  
474 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
475 *  
476 * @endcode
477 *
478 * Next is the typical loop over all cells to compute local contributions
479 * and then to transfer them into the global matrix and vector. The only
480 * change in this part, compared to @ref step_4 "step-4", is that we will use the
481 * <code>coefficient()</code> function defined above to compute the
482 * coefficient value at each quadrature point.
483 *
484 * @code
485 *   for (const auto &cell : dof_handler.active_cell_iterators())
486 *   {
487 *   fe_values.reinit(cell);
488 *  
489 *   cell_matrix = 0.;
490 *   cell_rhs = 0.;
491 *  
492 *   for (const unsigned int q_index : fe_values.quadrature_point_indices())
493 *   {
494 *   const double current_coefficient =
495 *   coefficient(fe_values.quadrature_point(q_index));
496 *   for (const unsigned int i : fe_values.dof_indices())
497 *   {
498 *   for (const unsigned int j : fe_values.dof_indices())
499 *   cell_matrix(i, j) +=
500 *   (current_coefficient * // a(x_q)
501 *   fe_values.shape_grad(i, q_index) * // grad phi_i(x_q)
502 *   fe_values.shape_grad(j, q_index) * // grad phi_j(x_q)
503 *   fe_values.JxW(q_index)); // dx
504 *  
505 *   cell_rhs(i) += (fe_values.shape_value(i, q_index) * // phi_i(x_q)
506 *   1.0 * // f(x_q)
507 *   fe_values.JxW(q_index)); // dx
508 *   }
509 *   }
510 *  
511 *  
512 *   cell->get_dof_indices(local_dof_indices);
513 *   for (const unsigned int i : fe_values.dof_indices())
514 *   {
515 *   for (const unsigned int j : fe_values.dof_indices())
516 *   system_matrix.add(local_dof_indices[i],
517 *   local_dof_indices[j],
518 *   cell_matrix(i, j));
519 *  
520 *   system_rhs(local_dof_indices[i]) += cell_rhs(i);
521 *   }
522 *   }
523 *  
524 * @endcode
525 *
526 * With the matrix so built, we use zero boundary values again:
527 *
528 * @code
529 *   std::map<types::global_dof_index, double> boundary_values;
531 *   types::boundary_id(0),
533 *   boundary_values);
534 *   MatrixTools::apply_boundary_values(boundary_values,
535 *   system_matrix,
536 *   solution,
537 *   system_rhs);
538 *   }
539 *  
540 *  
541 * @endcode
542 *
543 *
544 * <a name="step_5-Step5solve"></a>
545 * <h4>Step5::solve</h4>
546 *
547
548 *
549 * The solution process again looks mostly like in the previous
550 * examples. However, we will now use a preconditioned conjugate gradient
551 * algorithm. It is not very difficult to make this change. In fact, the only
552 * thing we have to alter is that we need an object which will act as a
553 * preconditioner. We will use SSOR (symmetric successive overrelaxation),
554 * with a relaxation factor of 1.2. For this purpose, the
555 * <code>SparseMatrix</code> class has a function which does one SSOR step,
556 * and we need to package the address of this function together with the
557 * matrix on which it should act (which is the matrix to be inverted) and the
558 * relaxation factor into one object. The <code>PreconditionSSOR</code> class
559 * does this for us. (<code>PreconditionSSOR</code> class takes a template
560 * argument denoting the matrix type it is supposed to work on. The default
561 * value is <code>SparseMatrix@<double@></code>, which is exactly what we need
562 * here, so we simply stick with the default and do not specify anything in
563 * the angle brackets.)
564 *
565
566 *
567 * Note that for the present case, SSOR doesn't really perform much better
568 * than most other preconditioners (though better than no preconditioning at
569 * all). A brief comparison of different preconditioners is presented in the
570 * Results section of the next tutorial program, @ref step_6 "step-6".
571 *
572
573 *
574 * With this, the rest of the function is trivial: instead of the
575 * <code>PreconditionIdentity</code> object we have created before, we now use
576 * the preconditioner we have declared, and the CG solver will do the rest for
577 * us:
578 *
579 * @code
580 *   template <int dim>
581 *   void Step5<dim>::solve()
582 *   {
583 *   SolverControl solver_control(1000, 1e-6 * system_rhs.l2_norm());
584 *   SolverCG<Vector<double>> solver(solver_control);
585 *  
586 *   PreconditionSSOR<SparseMatrix<double>> preconditioner;
587 *   preconditioner.initialize(system_matrix, 1.2);
588 *  
589 *   solver.solve(system_matrix, solution, system_rhs, preconditioner);
590 *  
591 *   std::cout << " " << solver_control.last_step()
592 *   << " CG iterations needed to obtain convergence." << std::endl;
593 *   }
594 *  
595 *  
596 * @endcode
597 *
598 *
599 * <a name="step_5-Step5output_resultsandsettingoutputflags"></a>
600 * <h4>Step5::output_results and setting output flags</h4>
601 *
602
603 *
604 * Writing output to a file is mostly the same as for the previous tutorial.
605 * The only difference is that we now need to construct a different filename
606 * for each refinement cycle.
607 *
608
609 *
610 * The function writes the output in VTU format, a variation of the VTK format
611 * that requires less disk space because it compresses the data. Of course,
612 * there are many other formats supported by the DataOut class if you
613 * desire to use a program for visualization that doesn't understand
614 * VTK or VTU.
615 *
616 * @code
617 *   template <int dim>
618 *   void Step5<dim>::output_results(const unsigned int cycle) const
619 *   {
620 *   DataOut<dim> data_out;
621 *  
622 *   data_out.attach_dof_handler(dof_handler);
623 *   data_out.add_data_vector(solution, "solution");
624 *  
625 *   data_out.build_patches();
626 *  
627 *   std::ofstream output("solution-" + std::to_string(cycle) + ".vtu");
628 *   data_out.write_vtu(output);
629 *   }
630 *  
631 *  
632 *  
633 * @endcode
634 *
635 *
636 * <a name="step_5-Step5run"></a>
637 * <h4>Step5::run</h4>
638 *
639
640 *
641 * The second to last thing in this program is the definition of the
642 * <code>run()</code> function. In contrast to the previous programs, we will
643 * compute on a sequence of meshes that after each iteration is globally
644 * refined. The function therefore consists of a loop over 6 cycles. In each
645 * cycle, we first print the cycle number, and then have to decide what to do
646 * with the mesh. If this is not the first cycle, we simply refine the
647 * existing mesh once globally. Before running through these cycles, however,
648 * we have to generate a mesh:
649 *
650
651 *
652 * In previous examples, we have already used some of the functions from the
653 * <code>GridGenerator</code> class. Here we would like to read a grid from a
654 * file where the cells are stored and which may originate from someone else,
655 * or may be the product of a mesh generator tool.
656 *
657
658 *
659 * In order to read a grid from a file, we generate an object of data type
660 * GridIn and associate the triangulation to it (i.e. we tell it to fill our
661 * triangulation object when we ask it to read the file). Then we open the
662 * respective file and initialize the triangulation with the data in the file :
663 *
664 * @code
665 *   template <int dim>
666 *   void Step5<dim>::run()
667 *   {
668 *   GridIn<dim> grid_in;
669 *   grid_in.attach_triangulation(triangulation);
670 *   std::ifstream input_file("circle-grid.inp");
671 * @endcode
672 *
673 * We would now like to read the file. However, the input file is only for a
674 * two-dimensional triangulation, while this function is a template for
675 * arbitrary dimension. Since this is only a demonstration program, we will
676 * not use different input files for the different dimensions, but rather
677 * quickly kill the whole program if we are not in 2d. Of course, since the
678 * main function below assumes that we are working in two dimensions we
679 * could skip this check, in this version of the program, without any ill
680 * effects.
681 *
682
683 *
684 * It turns out that perhaps 90 per cent of programming errors are invalid
685 * function parameters such as invalid array sizes, etc., so we use
686 * assertions heavily throughout deal.II to catch such mistakes. For this,
687 * the <code>Assert</code> macro is a good choice, since it makes sure that
688 * the condition which is given as first argument is valid, and if not
689 * throws an exception (its second argument) which will usually terminate
690 * the program giving information where the error occurred and what the
691 * reason was. (A longer discussion of what exactly the @p Assert macro
692 * does can be found in the @ref Exceptions "exception documentation topic".)
693 * This generally reduces the time to find programming errors
694 * dramatically and we have found assertions an invaluable means to program
695 * fast.
696 *
697
698 *
699 * On the other hand, all these checks (there are over 10,000 of them in the
700 * library at present) should not slow down the program too much if you want
701 * to do large computations. To this end, the <code>Assert</code> macro is
702 * only used in debug mode and expands to nothing if in optimized
703 * mode. Therefore, while you test your program on small problems and debug
704 * it, the assertions will tell you where the problems are. Once your
705 * program is stable, you can switch off debugging and the program will run
706 * your real computations without the assertions and at maximum speed. More
707 * precisely: turning off all the checks in the library (which prevent you
708 * from calling functions with wrong arguments, walking off of arrays, etc.)
709 * by compiling your program in optimized mode usually makes things run
710 * about four times faster. Even though optimized programs are more
711 * performant, you should always develop in debug mode since it allows
712 * the library to find lots of common programming errors automatically. For
713 * those who want to try: The way to switch from debug mode to optimized
714 * mode is to recompile your program with the command <code>make
715 * release</code>. The output of the <code>make</code> program should now
716 * indicate to you that the program is now compiled in optimized mode, and
717 * it will later also be linked to libraries that have been compiled for
718 * optimized mode. In order to switch back to debug mode, simply recompile
719 * with the command <code>make debug</code>.
720 *
721 * @code
722 *   Assert(dim == 2, ExcNotImplemented());
723 * @endcode
724 *
725 * ExcNotImplemented is a globally defined exception, which may be thrown
726 * whenever a piece of code has simply not been implemented for a case
727 * other than the condition checked in the assertion. Here, it would not
728 * be difficult to simply implement reading a *different* mesh file that
729 * contains a description of a 1d or 3d geometry, but this has not (yet)
730 * been implemented and so the exception is appropriate.
731 *
732
733 *
734 * Usually, one would like to use more specific
735 * exception classes, and particular in this case one would of course try
736 * to do something else if <code>dim</code> is not equal to two, e.g. create
737 * a grid using library functions. Aborting a program is usually not a good
738 * idea and assertions should really only be used for exceptional cases
739 * which should not occur, but might due to stupidity of the programmer,
740 * user, or someone else.
741 *
742
743 *
744 * So if we got past the assertion, we know that dim==2, and we can now
745 * actually read the grid. It is in UCD (unstructured cell data) format
746 * (though the convention is to use the suffix <code>inp</code> for UCD
747 * files):
748 *
749 * @code
750 *   grid_in.read_ucd(input_file);
751 * @endcode
752 *
753 * If you like to use another input format, you have to use one of the other
754 * <code>grid_in.read_xxx</code> function. (See the documentation of the
755 * <code>GridIn</code> class to find out what input formats are presently
756 * supported.)
757 *
758
759 *
760 * The grid in the file describes a circle. Therefore we have to use a
761 * manifold object which tells the triangulation where to put new points on
762 * the boundary when the grid is refined. Unlike @ref step_1 "step-1", since GridIn does
763 * not know that the domain has a circular boundary (unlike
764 * GridGenerator::hyper_shell) we have to explicitly attach a manifold to
765 * the boundary after creating the triangulation to get the correct result
766 * when we refine the mesh.
767 *
768 * @code
769 *   const SphericalManifold<dim> boundary;
770 *   triangulation.set_all_manifold_ids_on_boundary(0);
771 *   triangulation.set_manifold(0, boundary);
772 *  
773 *   for (unsigned int cycle = 0; cycle < 6; ++cycle)
774 *   {
775 *   std::cout << "Cycle " << cycle << ':' << std::endl;
776 *  
777 *   if (cycle != 0)
778 *   triangulation.refine_global(1);
779 *  
780 * @endcode
781 *
782 * Now that we have a mesh for sure, we write some output and do all the
783 * things that we have already seen in the previous examples.
784 *
785 * @code
786 *   std::cout << " Number of active cells: "
787 *   << triangulation.n_active_cells()
788 *   << std::endl
789 *   << " Total number of cells: "
790 *   << triangulation.n_cells()
791 *   << std::endl;
792 *  
793 *   setup_system();
794 *   assemble_system();
795 *   solve();
796 *   output_results(cycle);
797 *   }
798 *   }
799 *  
800 *  
801 * @endcode
802 *
803 *
804 * <a name="step_5-Thecodemaincodefunction"></a>
805 * <h3>The <code>main</code> function</h3>
806 *
807
808 *
809 * The main function looks mostly like the one in the previous example, so we
810 * won't comment on it further:
811 *
812 * @code
813 *   int main()
814 *   {
815 *   Step5<2> laplace_problem_2d;
816 *   laplace_problem_2d.run();
817 *   return 0;
818 *   }
819 * @endcode
820<a name="step_5-Results"></a><h1>Results</h1>
821
822
823
824Here is the console output:
825@code
826Cycle 0:
827 Number of active cells: 20
828 Total number of cells: 20
829 Number of degrees of freedom: 25
830 8 CG iterations needed to obtain convergence.
831Cycle 1:
832 Number of active cells: 80
833 Total number of cells: 100
834 Number of degrees of freedom: 89
835 12 CG iterations needed to obtain convergence.
836Cycle 2:
837 Number of active cells: 320
838 Total number of cells: 420
839 Number of degrees of freedom: 337
840 21 CG iterations needed to obtain convergence.
841Cycle 3:
842 Number of active cells: 1280
843 Total number of cells: 1700
844 Number of degrees of freedom: 1313
845 38 CG iterations needed to obtain convergence.
846Cycle 4:
847 Number of active cells: 5120
848 Total number of cells: 6820
849 Number of degrees of freedom: 5185
850 70 CG iterations needed to obtain convergence.
851Cycle 5:
852 Number of active cells: 20480
853 Total number of cells: 27300
854 Number of degrees of freedom: 20609
855 136 CG iterations needed to obtain convergence.
856@endcode
857
858
859
860In each cycle, the number of cells quadruples and the number of CG
861iterations roughly doubles.
862Also, in each cycle, the program writes one output graphic file in VTU
863format. They are depicted in the following:
864
865<table width="100%">
866 <tr>
867 <td>
868 <img src="https://www.dealii.org/images/steps/developer/step-5.solution-0-r9.2.png" alt="">
869 </td>
870 <td>
871 <img src="https://www.dealii.org/images/steps/developer/step-5.solution-1-r9.2.png" alt="">
872 </td>
873 </tr>
874 <tr>
875 <td>
876 <img src="https://www.dealii.org/images/steps/developer/step-5.solution-2-r9.2.png" alt="">
877 </td>
878 <td>
879 <img src="https://www.dealii.org/images/steps/developer/step-5.solution-3-r9.2.png" alt="">
880 </td>
881 </tr>
882 <tr>
883 <td>
884 <img src="https://www.dealii.org/images/steps/developer/step-5.solution-4-r9.2.png" alt="">
885 </td>
886 <td>
887 <img src="https://www.dealii.org/images/steps/developer/step-5.solution-5-r9.2.png" alt="">
888 </td>
889 </tr>
890</table>
891
892
893
894Due to the variable coefficient (the curvature there is reduced by the
895same factor by which the coefficient is increased), the top region of
896the solution is flattened. The gradient of the solution is
897discontinuous along the interface, although this is not very clearly
898visible in the pictures above. We will look at this in more detail in
899the next example.
900
901The pictures also show that the solution computed by this program is
902actually pretty wrong on a very coarse mesh (its magnitude is
903wrong). That's because no numerical method guarantees that the solution
904on a coarse mesh is particularly accurate -- but we know that the
905solution <i>converges</i> to the exact solution, and indeed you can
906see how the solutions from one mesh to the next seem to not change
907very much any more at the end.
908 *
909 *
910<a name="step_5-PlainProg"></a>
911<h1> The plain program</h1>
912@include "step-5.cc"
913*/
void initialize(const MatrixType &A, const AdditionalData &parameters=AdditionalData())
Point< 2 > first
Definition grid_out.cc:4623
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.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
std::vector< index_type > data
Definition mpi.cc:735
@ 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
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 reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation