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
step-72.h
Go to the documentation of this file.
1) const
490 *   {
491 *   return std::sin(2 * numbers::PI * (p[0] + p[1]));
492 *   }
493 *  
494 *  
495 * @endcode
496 *
497 *
498 * <a name="step_72-ThecodeMinimalSurfaceProblemcodeclassimplementation"></a>
499 * <h3>The <code>MinimalSurfaceProblem</code> class implementation</h3>
500 *
501
502 *
503 *
504 * <a name="step_72-MinimalSurfaceProblemMinimalSurfaceProblem"></a>
505 * <h4>MinimalSurfaceProblem::MinimalSurfaceProblem</h4>
506 *
507
508 *
509 * There have been no changes made to the class constructor.
510 *
511 * @code
512 *   template <int dim>
513 *   MinimalSurfaceProblem<dim>::MinimalSurfaceProblem()
514 *   : dof_handler(triangulation)
515 *   , fe(2)
516 *   , quadrature_formula(fe.degree + 1)
517 *   {}
518 *  
519 *  
520 * @endcode
521 *
522 *
523 * <a name="step_72-MinimalSurfaceProblemsetup_system"></a>
524 * <h4>MinimalSurfaceProblem::setup_system</h4>
525 *
526
527 *
528 * There have been no changes made to the function that sets up the class
529 * data structures, namely the DoFHandler, the hanging node constraints
530 * applied to the problem, and the linear system.
531 *
532 * @code
533 *   template <int dim>
534 *   void MinimalSurfaceProblem<dim>::setup_system()
535 *   {
536 *   dof_handler.distribute_dofs(fe);
537 *   current_solution.reinit(dof_handler.n_dofs());
538 *   newton_update.reinit(dof_handler.n_dofs());
539 *   system_rhs.reinit(dof_handler.n_dofs());
540 *  
541 *   zero_constraints.clear();
543 *   0,
545 *   zero_constraints);
546 *   DoFTools::make_hanging_node_constraints(dof_handler, zero_constraints);
547 *   zero_constraints.close();
548 *  
549 *   nonzero_constraints.clear();
551 *   0,
552 *   BoundaryValues<dim>(),
553 *   nonzero_constraints);
554 *   nonzero_constraints.close();
555 *  
556 *   DynamicSparsityPattern dsp(dof_handler.n_dofs());
557 *   DoFTools::make_sparsity_pattern(dof_handler, dsp, zero_constraints);
558 *  
559 *   sparsity_pattern.copy_from(dsp);
560 *   system_matrix.reinit(sparsity_pattern);
561 *   }
562 *  
563 * @endcode
564 *
565 *
566 * <a name="step_72-Assemblingthelinearsystem"></a>
567 * <h4>Assembling the linear system</h4>
568 *
569
570 *
571 *
572 * <a name="step_72-Manualassembly"></a>
573 * <h5>Manual assembly</h5>
574 *
575
576 *
577 * The assembly functions are the interesting contributions to this tutorial.
578 * The assemble_system_unassisted() method implements exactly the same
579 * assembly function as is detailed in @ref step_15 "step-15", but in this instance we
580 * use the MeshWorker::mesh_loop() function to multithread the assembly
581 * process. The reason for doing this is quite simple: When using
582 * automatic differentiation, we know that there is to be some additional
583 * computational overhead incurred. In order to mitigate this performance
584 * loss, we'd like to take advantage of as many (easily available)
585 * computational resources as possible. The MeshWorker::mesh_loop() concept
586 * makes this a relatively straightforward task. At the same time, for the
587 * purposes of fair comparison, we need to do the same to the implementation
588 * that uses no assistance when computing the residual or its linearization.
589 * (The MeshWorker::mesh_loop() function is first discussed in @ref step_12 "step-12" and
590 * @ref step_16 "step-16", if you'd like to read up on it.)
591 *
592
593 *
594 * The steps required to implement the multithreading are the same between the
595 * three functions, so we'll use the assemble_system_unassisted() function
596 * as an opportunity to focus on the multithreading itself.
597 *
598 * @code
599 *   template <int dim>
600 *   void MinimalSurfaceProblem<dim>::assemble_system_unassisted()
601 *   {
602 *   system_matrix = 0;
603 *   system_rhs = 0;
604 *  
605 *   const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
606 *  
607 * @endcode
608 *
609 * The MeshWorker::mesh_loop() expects that we provide two exemplar data
610 * structures. The first, `ScratchData`, is to store all large data that
611 * is to be reused between threads. The `CopyData` will hold the
612 * contributions to the linear system that come from each cell. These
613 * independent matrix-vector pairs must be accumulated into the
614 * global linear system sequentially. Since we don't need anything
615 * on top of what the MeshWorker::ScratchData and MeshWorker::CopyData
616 * classes already provide, we use these exact class definitions for
617 * our problem. Note that we only require a single instance of a local
618 * matrix, local right-hand side vector, and cell degree of freedom index
619 * vector -- the MeshWorker::CopyData therefore has `1` for all three
620 * of its template arguments.
621 *
622 * @code
623 *   using ScratchData = MeshWorker::ScratchData<dim>;
624 *   using CopyData = MeshWorker::CopyData<1, 1, 1>;
625 *  
626 * @endcode
627 *
628 * We also need to know what type of iterator we'll be working with
629 * during assembly. For simplicity, we just ask the compiler to work
630 * this out for us using the decltype() specifier, knowing that we'll
631 * be iterating over active cells owned by the @p dof_handler.
632 *
633 * @code
634 *   using CellIteratorType = decltype(dof_handler.begin_active());
635 *  
636 * @endcode
637 *
638 * Here we initialize the exemplar data structures. Since we know that
639 * we need to compute the shape function gradients, weighted Jacobian,
640 * and the position of the quadrate points in real space, we pass these
641 * flags into the class constructor.
642 *
643 * @code
644 *   const ScratchData sample_scratch_data(fe,
645 *   quadrature_formula,
646 *   update_gradients |
648 *   update_JxW_values);
649 *   const CopyData sample_copy_data(dofs_per_cell);
650 *  
651 * @endcode
652 *
653 * Now we define a lambda function that will perform the assembly on
654 * a single cell. The three arguments are those that will be expected by
655 * MeshWorker::mesh_loop(), due to the arguments that we'll pass to that
656 * final call. We also capture the @p this pointer, which means that we'll
657 * have access to "this" (i.e., the current `MinimalSurfaceProblem<dim>`)
658 * class instance, and its private member data (since the lambda function is
659 * defined within a MinimalSurfaceProblem<dim> method).
660 *
661
662 *
663 * At the top of the function, we initialize the data structures
664 * that are dependent on the cell for which the work is being
665 * performed. Observe that the reinitialization call actually
666 * returns an instance to an FEValues object that is initialized
667 * and stored within (and, therefore, reused by) the
668 * `scratch_data` object.
669 *
670
671 *
672 * Similarly, we get aliases to the local matrix, local RHS
673 * vector, and local cell DoF indices from the `copy_data`
674 * instance that MeshWorker::mesh_loop() provides. We then
675 * initialize the cell DoF indices, knowing that the local matrix
676 * and vector are already correctly sized.
677 *
678 * @code
679 *   const auto cell_worker = [this](const CellIteratorType &cell,
680 *   ScratchData &scratch_data,
681 *   CopyData &copy_data) {
682 *   const auto &fe_values = scratch_data.reinit(cell);
683 *  
684 *   FullMatrix<double> &cell_matrix = copy_data.matrices[0];
685 *   Vector<double> &cell_rhs = copy_data.vectors[0];
686 *   std::vector<types::global_dof_index> &local_dof_indices =
687 *   copy_data.local_dof_indices[0];
688 *   cell->get_dof_indices(local_dof_indices);
689 *  
690 * @endcode
691 *
692 * For Newton's method, we require the gradient of the solution at the
693 * point about which the problem is being linearized.
694 *
695
696 *
697 * Once we have that, we can perform assembly for this cell in
698 * the usual way. One minor difference to @ref step_15 "step-15" is that we've
699 * used the (rather convenient) range-based loops to iterate
700 * over all quadrature points and degrees-of-freedom.
701 *
702 * @code
703 *   std::vector<Tensor<1, dim>> old_solution_gradients(
704 *   fe_values.n_quadrature_points);
705 *   fe_values.get_function_gradients(current_solution,
706 *   old_solution_gradients);
707 *  
708 *   for (const unsigned int q : fe_values.quadrature_point_indices())
709 *   {
710 *   const double coeff =
711 *   1.0 / std::sqrt(1.0 + old_solution_gradients[q] *
712 *   old_solution_gradients[q]);
713 *  
714 *   for (const unsigned int i : fe_values.dof_indices())
715 *   {
716 *   for (const unsigned int j : fe_values.dof_indices())
717 *   cell_matrix(i, j) +=
718 *   (((fe_values.shape_grad(i, q) // ((\nabla \phi_i
719 *   * coeff // * a_n
720 *   * fe_values.shape_grad(j, q)) // * \nabla \phi_j)
721 *   - // -
722 *   (fe_values.shape_grad(i, q) // (\nabla \phi_i
723 *   * coeff * coeff * coeff // * a_n^3
724 *   * (fe_values.shape_grad(j, q) // * (\nabla \phi_j
725 *   * old_solution_gradients[q]) // * \nabla u_n)
726 *   * old_solution_gradients[q])) // * \nabla u_n)))
727 *   * fe_values.JxW(q)); // * dx
728 *  
729 *   cell_rhs(i) -= (fe_values.shape_grad(i, q) // \nabla \phi_i
730 *   * coeff // * a_n
731 *   * old_solution_gradients[q] // * \nabla u_n
732 *   * fe_values.JxW(q)); // * dx
733 *   }
734 *   }
735 *   };
736 *  
737 * @endcode
738 *
739 * The second lambda function that MeshWorker::mesh_loop() requires is
740 * one that performs the task of accumulating the local contributions
741 * in the global linear system. That is precisely what this one does,
742 * and the details of the implementation have been seen before. The
743 * primary point to recognize is that the local contributions are stored
744 * in the `copy_data` instance that is passed into this function. This
745 * `copy_data` has been filled with data during @a some call to the
746 * `cell_worker`.
747 *
748 * @code
749 *   const auto copier = [this](const CopyData &copy_data) {
750 *   const FullMatrix<double> &cell_matrix = copy_data.matrices[0];
751 *   const Vector<double> &cell_rhs = copy_data.vectors[0];
752 *   const std::vector<types::global_dof_index> &local_dof_indices =
753 *   copy_data.local_dof_indices[0];
754 *  
755 *   zero_constraints.distribute_local_to_global(
756 *   cell_matrix, cell_rhs, local_dof_indices, system_matrix, system_rhs);
757 *   };
758 *  
759 * @endcode
760 *
761 * We have all of the required functions definitions in place, so
762 * now we call the MeshWorker::mesh_loop() to perform the actual
763 * assembly. We pass a flag as the last parameter which states
764 * that we only want to perform the assembly on the
765 * cells. Internally, MeshWorker::mesh_loop() then distributes the
766 * available work to different threads, making efficient use of
767 * the multiple cores almost all of today's processors have to
768 * offer.
769 *
770 * @code
771 *   MeshWorker::mesh_loop(dof_handler.active_cell_iterators(),
772 *   cell_worker,
773 *   copier,
774 *   sample_scratch_data,
775 *   sample_copy_data,
776 *   MeshWorker::assemble_own_cells);
777 *   }
778 *  
779 * @endcode
780 *
781 *
782 * <a name="step_72-Assemblyviadifferentiationoftheresidualvector"></a>
783 * <h5>Assembly via differentiation of the residual vector</h5>
784 *
785
786 *
787 * As outlined in the introduction, what we need to do for this
788 * second approach is implement the local contributions @f$F(U)^K@f$
789 * from cell @f$K@f$ to the residual vector, and then let the
790 * AD machinery deal with how to compute the
791 * derivatives @f$J(U)_{ij}^K=\frac{\partial F(U)^K_i}{\partial U_j}@f$
792 * from it.
793 *
794
795 *
796 * For the following, recall that
797 * @f[
798 * F(U)_i^K \dealcoloneq
799 * \int\limits_K\nabla \varphi_i \cdot \left[ \frac{1}{\sqrt{1+|\nabla
800 * u|^{2}}} \nabla u \right] \, dV ,
801 * @f]
802 * where @f$u(\mathbf x)=\sum_j U_j \varphi_j(\mathbf x)@f$.
803 *
804
805 *
806 * Let us see how this is implemented in practice:
807 *
808 * @code
809 *   template <int dim>
810 *   void MinimalSurfaceProblem<dim>::assemble_system_with_residual_linearization()
811 *   {
812 *   system_matrix = 0;
813 *   system_rhs = 0;
814 *  
815 *   const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
816 *  
817 *   using ScratchData = MeshWorker::ScratchData<dim>;
818 *   using CopyData = MeshWorker::CopyData<1, 1, 1>;
819 *   using CellIteratorType = decltype(dof_handler.begin_active());
820 *  
821 *   const ScratchData sample_scratch_data(fe,
822 *   quadrature_formula,
823 *   update_gradients |
825 *   update_JxW_values);
826 *   const CopyData sample_copy_data(dofs_per_cell);
827 *  
828 * @endcode
829 *
830 * We'll define up front the AD data structures that we'll be using,
831 * utilizing the techniques shown in @ref step_71 "step-71".
832 * In this case, we choose the helper class that will automatically compute
833 * the linearization of the finite element residual using Sacado forward
834 * automatic differentiation types. These number types can be used to
835 * compute first derivatives only. This is exactly what we want, because we
836 * know that we'll only be linearizing the residual, which means that we
837 * only need to compute first-order derivatives. The return values from the
838 * calculations are to be of type `double`.
839 *
840
841 *
842 * We also need an extractor to retrieve some data related to the field
843 * solution to the problem.
844 *
845 * @code
846 *   using ADHelper = Differentiation::AD::ResidualLinearization<
847 *   Differentiation::AD::NumberTypes::sacado_dfad,
848 *   double>;
849 *   using ADNumberType = typename ADHelper::ad_type;
850 *  
851 *   const FEValuesExtractors::Scalar u_fe(0);
852 *  
853 * @endcode
854 *
855 * With this, let us define the lambda function that will be used
856 * to compute the cell contributions to the Jacobian matrix and
857 * the right hand side:
858 *
859 * @code
860 *   const auto cell_worker = [&u_fe, this](const CellIteratorType &cell,
861 *   ScratchData &scratch_data,
862 *   CopyData &copy_data) {
863 *   const auto &fe_values = scratch_data.reinit(cell);
864 *   const unsigned int dofs_per_cell = fe_values.get_fe().n_dofs_per_cell();
865 *  
866 *   FullMatrix<double> &cell_matrix = copy_data.matrices[0];
867 *   Vector<double> &cell_rhs = copy_data.vectors[0];
868 *   std::vector<types::global_dof_index> &local_dof_indices =
869 *   copy_data.local_dof_indices[0];
870 *   cell->get_dof_indices(local_dof_indices);
871 *  
872 * @endcode
873 *
874 * We'll now create and initialize an instance of the AD helper class.
875 * To do this, we need to specify how many independent variables and
876 * dependent variables there are. The independent variables will be the
877 * number of local degrees of freedom that our solution vector has,
878 * i.e., the number @f$j@f$ in the per-element representation of the
879 * discretized solution vector
880 * @f$u (\mathbf{x})|_K = \sum\limits_{j} U^K_i \varphi_j(\mathbf{x})@f$
881 * that indicates how many solution coefficients are associated with
882 * each finite element. In deal.II, this equals
883 * FiniteElement::dofs_per_cell. The number of dependent variables will be
884 * the number of entries in the local residual vector that we will be
885 * forming. In this particular problem (like many others that employ the
886 * [standard Galerkin
887 * method](https://en.wikipedia.org/wiki/Galerkin_method)) the number of
888 * local solution coefficients matches the number of local residual
889 * equations.
890 *
891 * @code
892 *   const unsigned int n_independent_variables = local_dof_indices.size();
893 *   const unsigned int n_dependent_variables = dofs_per_cell;
894 *   ADHelper ad_helper(n_independent_variables, n_dependent_variables);
895 *  
896 * @endcode
897 *
898 * Next we inform the helper of the values of the solution, i.e., the
899 * actual values for @f$U_j@f$ about which we
900 * wish to linearize. As this is done on each element individually, we
901 * have to extract the solution coefficients from the global solution
902 * vector. In other words, we define all of those coefficients @f$U_j@f$
903 * where @f$j@f$ is a local degree of freedom as the independent variables
904 * that enter the computation of the vector @f$F(U)^{K}@f$ (the dependent
905 * function).
906 *
907
908 *
909 * Then we get the complete set of degree of freedom values as
910 * represented by auto-differentiable numbers. The operations
911 * performed with these variables are tracked by the AD library
912 * from this point until the object goes out of scope. So it is
913 * <em>precisely these variables</em> with respect to which we will
914 * compute derivatives of the residual entries.
915 *
916 * @code
917 *   ad_helper.register_dof_values(current_solution, local_dof_indices);
918 *  
919 *   const std::vector<ADNumberType> &dof_values_ad =
920 *   ad_helper.get_sensitive_dof_values();
921 *  
922 * @endcode
923 *
924 * Then we do some problem specific tasks, the first being to
925 * compute all values, (spatial) gradients, and the like based on
926 * "sensitive" AD degree of freedom values. In this instance we are
927 * retrieving the solution gradients at each quadrature point. Observe
928 * that the solution gradients are now sensitive
929 * to the values of the degrees of freedom as they use the @p ADNumberType
930 * as the scalar type and the @p dof_values_ad vector provides the local
931 * DoF values.
932 *
933 * @code
934 *   std::vector<Tensor<1, dim, ADNumberType>> old_solution_gradients(
935 *   fe_values.n_quadrature_points);
936 *   fe_values[u_fe].get_function_gradients_from_local_dof_values(
937 *   dof_values_ad, old_solution_gradients);
938 *  
939 * @endcode
940 *
941 * The next variable that we declare will store the cell residual vector
942 * contributions. This is rather self-explanatory, save for one
943 * <b>very important</b> detail:
944 * Note that each entry in the vector is hand-initialized with a value
945 * of zero. This is a <em>highly recommended</em> practice, as some AD
946 * libraries appear not to safely initialize the internal data
947 * structures of these number types. Not doing so could lead to some
948 * very hard to understand or detect bugs (appreciate that the author
949 * of this program mentions this out of, generally bad, experience). So
950 * out of an abundance of caution it's worthwhile zeroing the initial
951 * value explicitly. After that, apart from a sign change the residual
952 * assembly looks much the same as we saw for the cell RHS vector before:
953 * We loop over all quadrature points, ensure that the coefficient now
954 * encodes its dependence on the (sensitive) finite element DoF values by
955 * using the correct `ADNumberType`, and finally we assemble the
956 * components of the residual vector. For complete clarity, the finite
957 * element shape functions (and their gradients, etc.) as well as the
958 * "JxW" values remain scalar
959 * valued, but the @p coeff and the @p old_solution_gradients at each
960 * quadrature point are computed in terms of the independent
961 * variables.
962 *
963 * @code
964 *   std::vector<ADNumberType> residual_ad(n_dependent_variables,
965 *   ADNumberType(0.0));
966 *   for (const unsigned int q : fe_values.quadrature_point_indices())
967 *   {
968 *   const ADNumberType coeff =
969 *   1.0 / std::sqrt(1.0 + old_solution_gradients[q] *
970 *   old_solution_gradients[q]);
971 *  
972 *   for (const unsigned int i : fe_values.dof_indices())
973 *   {
974 *   residual_ad[i] += (fe_values.shape_grad(i, q) // \nabla \phi_i
975 *   * coeff // * a_n
976 *   * old_solution_gradients[q]) // * \nabla u_n
977 *   * fe_values.JxW(q); // * dx
978 *   }
979 *   }
980 *  
981 * @endcode
982 *
983 * Once we have the full cell residual vector computed, we can register
984 * it with the helper class.
985 *
986
987 *
988 * Thereafter, we compute the residual values (basically,
989 * extracting the real values from what we already computed) and
990 * their Jacobian (the linearization of each residual component
991 * with respect to all cell DoFs) at the evaluation point. For
992 * the purposes of assembly into the global linear system, we
993 * have to respect the sign difference between the residual and
994 * the RHS contribution: For Newton's method, the right hand
995 * side vector needs to be equal to the *negative* residual
996 * vector.
997 *
998 * @code
999 *   ad_helper.register_residual_vector(residual_ad);
1000 *  
1001 *   ad_helper.compute_residual(cell_rhs);
1002 *   cell_rhs *= -1.0;
1003 *  
1004 *   ad_helper.compute_linearization(cell_matrix);
1005 *   };
1006 *  
1007 * @endcode
1008 *
1009 * The remainder of the function equals what we had previously:
1010 *
1011 * @code
1012 *   const auto copier = [this](const CopyData &copy_data) {
1013 *   const FullMatrix<double> &cell_matrix = copy_data.matrices[0];
1014 *   const Vector<double> &cell_rhs = copy_data.vectors[0];
1015 *   const std::vector<types::global_dof_index> &local_dof_indices =
1016 *   copy_data.local_dof_indices[0];
1017 *  
1018 *   zero_constraints.distribute_local_to_global(
1019 *   cell_matrix, cell_rhs, local_dof_indices, system_matrix, system_rhs);
1020 *   };
1021 *  
1022 *   MeshWorker::mesh_loop(dof_handler.active_cell_iterators(),
1023 *   cell_worker,
1024 *   copier,
1025 *   sample_scratch_data,
1026 *   sample_copy_data,
1028 *   }
1029 *  
1030 * @endcode
1031 *
1032 *
1033 * <a name="step_72-Assemblyviadifferentiationoftheenergyfunctional"></a>
1034 * <h5>Assembly via differentiation of the energy functional</h5>
1035 *
1036
1037 *
1038 * In this third approach, we compute residual and Jacobian as first
1039 * and second derivatives of the local energy functional
1040 * @f[
1041 * E\left( U \right)^K
1042 * \dealcoloneq \int\limits_{K} \Psi \left( u \right) \, dV
1043 * \approx \sum\limits_{q}^{n_{\textrm{q-points}}} \Psi \left( u \left(
1044 * \mathbf{X}_{q} \right) \right) \underbrace{\vert J_{q} \vert \times
1045 * W_{q}}_{\text{JxW(q)}}
1046 * @f]
1047 * with the energy density given by
1048 * @f[
1049 * \Psi \left( u \right) = \sqrt{1+|\nabla u|^{2}} .
1050 * @f]
1051 *
1052
1053 *
1054 * Let us again see how this is done:
1055 *
1056 * @code
1057 *   template <int dim>
1058 *   void MinimalSurfaceProblem<dim>::assemble_system_using_energy_functional()
1059 *   {
1060 *   system_matrix = 0;
1061 *   system_rhs = 0;
1062 *  
1063 *   const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
1064 *  
1065 *   using ScratchData = MeshWorker::ScratchData<dim>;
1066 *   using CopyData = MeshWorker::CopyData<1, 1, 1>;
1067 *   using CellIteratorType = decltype(dof_handler.begin_active());
1068 *  
1069 *   const ScratchData sample_scratch_data(fe,
1070 *   quadrature_formula,
1071 *   update_gradients |
1073 *   update_JxW_values);
1074 *   const CopyData sample_copy_data(dofs_per_cell);
1075 *  
1076 * @endcode
1077 *
1078 * In this implementation of the assembly process, we choose the helper
1079 * class that will automatically compute both the residual and its
1080 * linearization from the cell contribution to an energy functional using
1081 * nested Sacado forward automatic differentiation types.
1082 * The selected number types can be used to compute both first and
1083 * second derivatives. We require this, as the residual defined as the
1084 * sensitivity of the potential energy with respect to the DoF values (i.e.
1085 * its gradient). We'll then need to linearize the residual, implying that
1086 * second derivatives of the potential energy must be computed. You might
1087 * want to compare this with the definition of `ADHelper` used int
1088 * previous function, where we used
1089 * `Differentiation::AD::ResidualLinearization<Differentiation::AD::NumberTypes::sacado_dfad,double>`.
1090 *
1091 * @code
1092 *   using ADHelper = Differentiation::AD::EnergyFunctional<
1093 *   Differentiation::AD::NumberTypes::sacado_dfad_dfad,
1094 *   double>;
1095 *   using ADNumberType = typename ADHelper::ad_type;
1096 *  
1097 *   const FEValuesExtractors::Scalar u_fe(0);
1098 *  
1099 * @endcode
1100 *
1101 * Let us then again define the lambda function that does the integration on
1102 * a cell.
1103 *
1104
1105 *
1106 * To initialize an instance of the helper class, we now only require
1107 * that the number of independent variables (that is, the number
1108 * of degrees of freedom associated with the element solution vector)
1109 * are known up front. This is because the second-derivative matrix that
1110 * results from an energy functional is necessarily square (and also,
1111 * incidentally, symmetric).
1112 *
1113 * @code
1114 *   const auto cell_worker = [&u_fe, this](const CellIteratorType &cell,
1115 *   ScratchData &scratch_data,
1116 *   CopyData &copy_data) {
1117 *   const auto &fe_values = scratch_data.reinit(cell);
1118 *  
1119 *   FullMatrix<double> &cell_matrix = copy_data.matrices[0];
1120 *   Vector<double> &cell_rhs = copy_data.vectors[0];
1121 *   std::vector<types::global_dof_index> &local_dof_indices =
1122 *   copy_data.local_dof_indices[0];
1123 *   cell->get_dof_indices(local_dof_indices);
1124 *  
1125 *   const unsigned int n_independent_variables = local_dof_indices.size();
1126 *   ADHelper ad_helper(n_independent_variables);
1127 *  
1128 * @endcode
1129 *
1130 * Once more, we register all cell DoFs values with the helper, followed
1131 * by extracting the "sensitive" variant of these values that are to be
1132 * used in subsequent operations that must be differentiated -- one of
1133 * those being the calculation of the solution gradients.
1134 *
1135 * @code
1136 *   ad_helper.register_dof_values(current_solution, local_dof_indices);
1137 *  
1138 *   const std::vector<ADNumberType> &dof_values_ad =
1139 *   ad_helper.get_sensitive_dof_values();
1140 *  
1141 *   std::vector<Tensor<1, dim, ADNumberType>> old_solution_gradients(
1142 *   fe_values.n_quadrature_points);
1143 *   fe_values[u_fe].get_function_gradients_from_local_dof_values(
1144 *   dof_values_ad, old_solution_gradients);
1145 *  
1146 * @endcode
1147 *
1148 * We next create a variable that stores the cell total energy.
1149 * Once more we emphasize that we explicitly zero-initialize this value,
1150 * thereby ensuring the integrity of the data for this starting value.
1151 *
1152
1153 *
1154 * The aim for our approach is then to compute the cell total
1155 * energy, which is the sum of the internal (due to right hand
1156 * side functions, typically linear in @f$U@f$) and external
1157 * energies. In this particular case, we have no external
1158 * energies (e.g., from source terms or Neumann boundary
1159 * conditions), so we'll focus on the internal energy part.
1160 *
1161
1162 *
1163 * In fact, computing @f$E(U)^K@f$ is almost trivial, requiring only
1164 * the following lines:
1165 *
1166 * @code
1167 *   ADNumberType energy_ad = ADNumberType(0.0);
1168 *   for (const unsigned int q : fe_values.quadrature_point_indices())
1169 *   {
1170 *   const ADNumberType psi = std::sqrt(1.0 + old_solution_gradients[q] *
1171 *   old_solution_gradients[q]);
1172 *  
1173 *   energy_ad += psi * fe_values.JxW(q);
1174 *   }
1175 *  
1176 * @endcode
1177 *
1178 * After we've computed the total energy on this cell, we'll
1179 * register it with the helper. Based on that, we may now
1180 * compute the desired quantities, namely the residual values
1181 * and their Jacobian at the evaluation point. As before, the
1182 * Newton right hand side needs to be the negative of the
1183 * residual:
1184 *
1185 * @code
1186 *   ad_helper.register_energy_functional(energy_ad);
1187 *  
1188 *   ad_helper.compute_residual(cell_rhs);
1189 *   cell_rhs *= -1.0;
1190 *  
1191 *   ad_helper.compute_linearization(cell_matrix);
1192 *   };
1193 *  
1194 * @endcode
1195 *
1196 * As in the previous two functions, the remainder of the function is as
1197 * before:
1198 *
1199 * @code
1200 *   const auto copier = [this](const CopyData &copy_data) {
1201 *   const FullMatrix<double> &cell_matrix = copy_data.matrices[0];
1202 *   const Vector<double> &cell_rhs = copy_data.vectors[0];
1203 *   const std::vector<types::global_dof_index> &local_dof_indices =
1204 *   copy_data.local_dof_indices[0];
1205 *  
1206 *   zero_constraints.distribute_local_to_global(
1207 *   cell_matrix, cell_rhs, local_dof_indices, system_matrix, system_rhs);
1208 *   };
1209 *  
1210 *   MeshWorker::mesh_loop(dof_handler.active_cell_iterators(),
1211 *   cell_worker,
1212 *   copier,
1213 *   sample_scratch_data,
1214 *   sample_copy_data,
1216 *   }
1217 *  
1218 *  
1219 * @endcode
1220 *
1221 *
1222 * <a name="step_72-MinimalSurfaceProblemsolve"></a>
1223 * <h4>MinimalSurfaceProblem::solve</h4>
1224 *
1225
1226 *
1227 * The solve function is the same as is used in @ref step_15 "step-15".
1228 *
1229 * @code
1230 *   template <int dim>
1231 *   void MinimalSurfaceProblem<dim>::solve()
1232 *   {
1233 *   SolverControl solver_control(system_rhs.size(),
1234 *   system_rhs.l2_norm() * 1e-6);
1235 *   SolverCG<Vector<double>> solver(solver_control);
1236 *  
1237 *   PreconditionSSOR<SparseMatrix<double>> preconditioner;
1238 *   preconditioner.initialize(system_matrix, 1.2);
1239 *  
1240 *   solver.solve(system_matrix, newton_update, system_rhs, preconditioner);
1241 *  
1242 *   zero_constraints.distribute(newton_update);
1243 *  
1244 *   const double alpha = determine_step_length();
1245 *   current_solution.add(alpha, newton_update);
1246 *   }
1247 *  
1248 *  
1249 * @endcode
1250 *
1251 *
1252 * <a name="step_72-MinimalSurfaceProblemrefine_mesh"></a>
1253 * <h4>MinimalSurfaceProblem::refine_mesh</h4>
1254 *
1255
1256 *
1257 * Nothing has changed since @ref step_15 "step-15" with respect to the mesh refinement
1258 * procedure and transfer of the solution between adapted meshes.
1259 *
1260 * @code
1261 *   template <int dim>
1262 *   void MinimalSurfaceProblem<dim>::refine_mesh()
1263 *   {
1264 *   Vector<float> estimated_error_per_cell(triangulation.n_active_cells());
1265 *  
1267 *   dof_handler,
1268 *   QGauss<dim - 1>(fe.degree + 1),
1269 *   std::map<types::boundary_id, const Function<dim> *>(),
1270 *   current_solution,
1271 *   estimated_error_per_cell);
1272 *  
1274 *   estimated_error_per_cell,
1275 *   0.3,
1276 *   0.03);
1277 *  
1278 *   triangulation.prepare_coarsening_and_refinement();
1279 *  
1280 *   SolutionTransfer<dim> solution_transfer(dof_handler);
1281 *   const Vector<double> coarse_solution = current_solution;
1282 *   solution_transfer.prepare_for_coarsening_and_refinement(coarse_solution);
1283 *  
1284 *   triangulation.execute_coarsening_and_refinement();
1285 *  
1286 *   setup_system();
1287 *  
1288 *   solution_transfer.interpolate(current_solution);
1289 *   nonzero_constraints.distribute(current_solution);
1290 *   }
1291 *  
1292 *  
1293 *  
1294 * @endcode
1295 *
1296 *
1297 * <a name="step_72-MinimalSurfaceProblemcompute_residual"></a>
1298 * <h4>MinimalSurfaceProblem::compute_residual</h4>
1299 *
1300
1301 *
1302 * ... as does the function used to compute the residual during the
1303 * solution iteration procedure. One could replace this by
1304 * differentiation of the energy functional if one really wanted,
1305 * but for simplicity we here simply copy what we already had in
1306 * @ref step_15 "step-15".
1307 *
1308 * @code
1309 *   template <int dim>
1310 *   double MinimalSurfaceProblem<dim>::compute_residual(const double alpha) const
1311 *   {
1312 *   Vector<double> residual(dof_handler.n_dofs());
1313 *  
1314 *   Vector<double> evaluation_point(dof_handler.n_dofs());
1315 *   evaluation_point = current_solution;
1316 *   evaluation_point.add(alpha, newton_update);
1317 *  
1318 *   FEValues<dim> fe_values(fe,
1319 *   quadrature_formula,
1321 *   update_JxW_values);
1322 *  
1323 *   const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
1324 *   const unsigned int n_q_points = quadrature_formula.size();
1325 *  
1326 *   Vector<double> cell_residual(dofs_per_cell);
1327 *   std::vector<Tensor<1, dim>> gradients(n_q_points);
1328 *  
1329 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1330 *  
1331 *   for (const auto &cell : dof_handler.active_cell_iterators())
1332 *   {
1333 *   cell_residual = 0;
1334 *   fe_values.reinit(cell);
1335 *  
1336 *   fe_values.get_function_gradients(evaluation_point, gradients);
1337 *  
1338 *   for (unsigned int q = 0; q < n_q_points; ++q)
1339 *   {
1340 *   const double coeff =
1341 *   1.0 / std::sqrt(1.0 + gradients[q] * gradients[q]);
1342 *  
1343 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1344 *   cell_residual(i) -= (fe_values.shape_grad(i, q) // \nabla \phi_i
1345 *   * coeff // * a_n
1346 *   * gradients[q] // * \nabla u_n
1347 *   * fe_values.JxW(q)); // * dx
1348 *   }
1349 *  
1350 *   cell->get_dof_indices(local_dof_indices);
1351 *   zero_constraints.distribute_local_to_global(cell_residual,
1352 *   local_dof_indices,
1353 *   residual);
1354 *   }
1355 *  
1356 *   return residual.l2_norm();
1357 *   }
1358 *  
1359 *  
1360 *  
1361 * @endcode
1362 *
1363 *
1364 * <a name="step_72-MinimalSurfaceProblemdetermine_step_length"></a>
1365 * <h4>MinimalSurfaceProblem::determine_step_length</h4>
1366 *
1367
1368 *
1369 * The choice of step length (or, under-relaxation factor) for the nonlinear
1370 * iterations procedure remains fixed at the value chosen and discussed in
1371 * @ref step_15 "step-15".
1372 *
1373 * @code
1374 *   template <int dim>
1375 *   double MinimalSurfaceProblem<dim>::determine_step_length() const
1376 *   {
1377 *   return 0.1;
1378 *   }
1379 *  
1380 *  
1381 *  
1382 * @endcode
1383 *
1384 *
1385 * <a name="step_72-MinimalSurfaceProblemoutput_results"></a>
1386 * <h4>MinimalSurfaceProblem::output_results</h4>
1387 *
1388
1389 *
1390 * This last function to be called from `run()` outputs the current solution
1391 * (and the Newton update) in graphical form as a VTU file. It is entirely the
1392 * same as what has been used in previous tutorials.
1393 *
1394 * @code
1395 *   template <int dim>
1396 *   void MinimalSurfaceProblem<dim>::output_results(
1397 *   const unsigned int refinement_cycle) const
1398 *   {
1399 *   DataOut<dim> data_out;
1400 *  
1401 *   data_out.attach_dof_handler(dof_handler);
1402 *   data_out.add_data_vector(current_solution, "solution");
1403 *   data_out.add_data_vector(newton_update, "update");
1404 *   data_out.build_patches();
1405 *  
1406 *   const std::string filename =
1407 *   "solution-" + Utilities::int_to_string(refinement_cycle, 2) + ".vtu";
1408 *   std::ofstream output(filename);
1409 *   data_out.write_vtu(output);
1410 *   }
1411 *  
1412 *  
1413 * @endcode
1414 *
1415 *
1416 * <a name="step_72-MinimalSurfaceProblemrun"></a>
1417 * <h4>MinimalSurfaceProblem::run</h4>
1418 *
1419
1420 *
1421 * In the run function, most remains the same as was first implemented
1422 * in @ref step_15 "step-15". The only observable changes are that we can now choose (via
1423 * the parameter file) what the final acceptable tolerance for the system
1424 * residual is, and that we can choose which method of assembly we wish to
1425 * utilize. To make the second choice clear, we output to the console some
1426 * message which indicates the selection. Since we're interested in comparing
1427 * the time taken to assemble for each of the three methods, we've also
1428 * added a timer that keeps a track of how much time is spent during assembly.
1429 * We also track the time taken to solve the linear system, so that we can
1430 * contrast those numbers to the part of the code which would normally take
1431 * the longest time to execute.
1432 *
1433 * @code
1434 *   template <int dim>
1435 *   void MinimalSurfaceProblem<dim>::run(const int formulation,
1436 *   const double tolerance)
1437 *   {
1438 *   std::cout << "******** Assembly approach ********" << std::endl;
1439 *   const std::array<std::string, 3> method_descriptions = {
1440 *   {"Unassisted implementation (full hand linearization).",
1441 *   "Automated linearization of the finite element residual.",
1442 *   "Automated computation of finite element residual and linearization using a variational formulation."}};
1443 *   AssertIndexRange(formulation, method_descriptions.size());
1444 *   std::cout << method_descriptions[formulation] << std::endl << std::endl;
1445 *  
1446 *  
1448 *  
1449 *   GridGenerator::hyper_ball(triangulation);
1450 *   triangulation.refine_global(2);
1451 *  
1452 *   setup_system();
1453 *   nonzero_constraints.distribute(current_solution);
1454 *  
1455 *   double last_residual_norm = std::numeric_limits<double>::max();
1456 *   unsigned int refinement_cycle = 0;
1457 *   do
1458 *   {
1459 *   std::cout << "Mesh refinement step " << refinement_cycle << std::endl;
1460 *  
1461 *   if (refinement_cycle != 0)
1462 *   refine_mesh();
1463 *  
1464 *   std::cout << " Initial residual: " << compute_residual(0) << std::endl;
1465 *  
1466 *   for (unsigned int inner_iteration = 0; inner_iteration < 5;
1467 *   ++inner_iteration)
1468 *   {
1469 *   {
1470 *   TimerOutput::Scope t(timer, "Assemble");
1471 *  
1472 *   if (formulation == 0)
1473 *   assemble_system_unassisted();
1474 *   else if (formulation == 1)
1475 *   assemble_system_with_residual_linearization();
1476 *   else if (formulation == 2)
1477 *   assemble_system_using_energy_functional();
1478 *   else
1479 *   AssertThrow(false, ExcNotImplemented());
1480 *   }
1481 *  
1482 *   last_residual_norm = system_rhs.l2_norm();
1483 *  
1484 *   {
1485 *   TimerOutput::Scope t(timer, "Solve");
1486 *   solve();
1487 *   }
1488 *  
1489 *  
1490 *   std::cout << " Residual: " << compute_residual(0) << std::endl;
1491 *   }
1492 *  
1493 *   output_results(refinement_cycle);
1494 *  
1495 *   ++refinement_cycle;
1496 *   std::cout << std::endl;
1497 *   }
1498 *   while (last_residual_norm > tolerance);
1499 *   }
1500 *   } // namespace Step72
1501 *  
1502 * @endcode
1503 *
1504 *
1505 * <a name="step_72-Themainfunction"></a>
1506 * <h4>The main function</h4>
1507 *
1508
1509 *
1510 * Finally the main function. This follows the scheme of most other main
1511 * functions, with two obvious exceptions:
1512 * - We call Utilities::MPI::MPI_InitFinalize in order to set up (via a hidden
1513 * default parameter) the number of threads using the execution of
1514 * multithreaded tasks.
1515 * - We also have a few lines dedicates to reading in or initializing the
1516 * user-defined parameters that will be considered during the execution of the
1517 * program.
1518 *
1519 * @code
1520 *   int main(int argc, char *argv[])
1521 *   {
1522 *   try
1523 *   {
1524 *   using namespace Step72;
1525 *  
1526 *   Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv);
1527 *  
1528 *   std::string prm_file;
1529 *   if (argc > 1)
1530 *   prm_file = argv[1];
1531 *   else
1532 *   prm_file = "parameters.prm";
1533 *  
1534 *   const MinimalSurfaceProblemParameters parameters;
1535 *   ParameterAcceptor::initialize(prm_file);
1536 *  
1537 *   MinimalSurfaceProblem<2> minimal_surface_problem_2d;
1538 *   minimal_surface_problem_2d.run(parameters.formulation,
1539 *   parameters.tolerance);
1540 *   }
1541 *   catch (std::exception &exc)
1542 *   {
1543 *   std::cerr << std::endl
1544 *   << std::endl
1545 *   << "----------------------------------------------------"
1546 *   << std::endl;
1547 *   std::cerr << "Exception on processing: " << std::endl
1548 *   << exc.what() << std::endl
1549 *   << "Aborting!" << std::endl
1550 *   << "----------------------------------------------------"
1551 *   << std::endl;
1552 *  
1553 *   return 1;
1554 *   }
1555 *   catch (...)
1556 *   {
1557 *   std::cerr << std::endl
1558 *   << std::endl
1559 *   << "----------------------------------------------------"
1560 *   << std::endl;
1561 *   std::cerr << "Unknown exception!" << std::endl
1562 *   << "Aborting!" << std::endl
1563 *   << "----------------------------------------------------"
1564 *   << std::endl;
1565 *   return 1;
1566 *   }
1567 *   return 0;
1568 *   }
1569 * @endcode
1570<a name="step_72-Results"></a><h1>Results</h1>
1571
1572
1573Since there was no change to the physics of the problem that has first been analyzed
1574in @ref step_15 "step-15", there is nothing to report about that. The only outwardly noticeable
1575difference between them is that, by default, this program will only run 9 mesh
1576refinement steps (as opposed to @ref step_15 "step-15", which executes 11 refinements).
1577This will be observable in the simulation status that appears between the
1578header text that prints which assembly method is being used, and the final
1579timings. (All timings reported below were obtained in release mode.)
1580
1581@code
1582Mesh refinement step 0
1583 Initial residual: 1.53143
1584 Residual: 1.08746
1585 Residual: 0.966748
1586 Residual: 0.859602
1587 Residual: 0.766462
1588 Residual: 0.685475
1589
1590...
1591
1592Mesh refinement step 9
1593 Initial residual: 0.00924594
1594 Residual: 0.00831928
1595 Residual: 0.0074859
1596 Residual: 0.0067363
1597 Residual: 0.00606197
1598 Residual: 0.00545529
1599@endcode
1600
1601So what is interesting for us to compare is how long the assembly process takes
1602for the three different implementations, and to put that into some greater context.
1603Below is the output for the hand linearization (as computed on a circa 2012
1604four core, eight thread laptop -- but we're only really interested in the
1605relative time between the different implementations):
1606@code
1607******** Assembly approach ********
1608Unassisted implementation (full hand linearization).
1609
1610...
1611
1612+---------------------------------------------+------------+------------+
1613| Total wallclock time elapsed since start | 35.1s | |
1614| | | |
1615| Section | no. calls | wall time | % of total |
1616+---------------------------------+-----------+------------+------------+
1617| Assemble | 50 | 1.56s | 4.5% |
1618| Solve | 50 | 30.8s | 88% |
1619+---------------------------------+-----------+------------+------------+
1620@endcode
1621And for the implementation that linearizes the residual in an automated
1622manner using the Sacado dynamic forward AD number type:
1623@code
1624******** Assembly approach ********
1625Automated linearization of the finite element residual.
1626
1627...
1628
1629+---------------------------------------------+------------+------------+
1630| Total wallclock time elapsed since start | 40.1s | |
1631| | | |
1632| Section | no. calls | wall time | % of total |
1633+---------------------------------+-----------+------------+------------+
1634| Assemble | 50 | 8.8s | 22% |
1635| Solve | 50 | 28.6s | 71% |
1636+---------------------------------+-----------+------------+------------+
1637@endcode
1638And, lastly, for the implementation that computes both the residual and
1639its linearization directly from an energy functional (using nested Sacado
1640dynamic forward AD numbers):
1641@code
1642******** Assembly approach ********
1643Automated computation of finite element residual and linearization using a variational formulation.
1644
1645...
1646
1647+---------------------------------------------+------------+------------+
1648| Total wallclock time elapsed since start | 48.8s | |
1649| | | |
1650| Section | no. calls | wall time | % of total |
1651+---------------------------------+-----------+------------+------------+
1652| Assemble | 50 | 16.7s | 34% |
1653| Solve | 50 | 29.3s | 60% |
1654+---------------------------------+-----------+------------+------------+
1655@endcode
1656
1657It's evident that the more work that is passed off to the automatic differentiation
1658framework to perform, the more time is spent during the assembly process. Accumulated
1659over all refinement steps, using one level of automatic differentiation resulted
1660in @f$5.65 \times@f$ more computational time being spent in the assembly stage when
1661compared to unassisted assembly, while assembling the discrete linear system took
1662@f$10.7 \times@f$ longer when deriving directly from the energy functional.
1663Unsurprisingly, the overall time spent solving the linear system remained unchanged.
1664This means that the proportion of time spent in the solve phase to the assembly phase
1665shifted significantly as the number of times automated differentiation was performed
1666at the finite element level. For many, this might mean that leveraging higher-order
1667differentiation (at the finite element level) in production code leads to an
1668unacceptable overhead, but it may still be useful during the prototyping phase.
1669A good compromise between the two may, therefore, be the automated linearization
1670of the finite element residual, which offers a lot of convenience at a measurable,
1671but perhaps not unacceptable, cost. Alternatively, one could consider
1672not re-building the Newton matrix in every step -- a topic that is
1673explored in substantial depth in @ref step_77 "step-77".
1674
1675Of course, in practice the actual overhead is very much dependent on the problem being evaluated
1676(e.g., how many components there are in the solution, what the nature of the function
1677being differentiated is, etc.). So the exact results presented here should be
1678interpreted within the context of this scalar problem alone, and when it comes to
1679other problems, some preliminary investigation by the user is certainly warranted.
1680
1681
1682<a name="step_72-Possibilitiesforextensions"></a><h3> Possibilities for extensions </h3>
1683
1684
1685Like @ref step_71 "step-71", there are a few items related to automatic differentiation that could
1686be evaluated further:
1687- The use of other AD frameworks should be investigated, with the outlook that
1688 alternative implementations may provide performance benefits.
1689- It is also worth evaluating AD number types other than those that have been
1690 hard-coded into this tutorial. With regard to twice differentiable types
1691 employed at the finite-element level, mixed differentiation modes ("RAD-FAD")
1692 should in principle be more computationally efficient than the single
1693 mode ("FAD-FAD") types employed here. The reason that the RAD-FAD type was not
1694 selected by default is that, at the time of writing, there remain some
1695 bugs in its implementation within the Sacado library that lead to memory leaks.
1696 This is documented in the @ref auto_symb_diff topic.
1697- It might be the case that using reduced precision types (i.e., `float`) as the
1698 scalar types for the AD numbers could render a reduction in computational
1699 expense during assembly. Using `float` as the data type for the
1700 matrix and the residual is not unreasonable, given that the Newton
1701 update is only meant to get us closer to the solution, but not
1702 actually *to* the solution; as a consequence, it makes sense to
1703 consider using reduced-precision data types for computing these
1704 updates, and then accumulating these updates in a solution vector
1705 that uses the full `double` precision accuracy.
1706- One further method of possibly reducing resources during assembly is to frame
1707 the AD implementations as a constitutive model. This would be similar to the
1708 approach adopted in @ref step_71 "step-71", and pushes the starting point for the automatic
1709 differentiation one level higher up the chain of computations. This, in turn,
1710 means that less operations are tracked by the AD library, thereby reducing the
1711 cost of differentiating (even though one would perform the differentiation at
1712 each cell quadrature point).
1713- @ref step_77 "step-77" is yet another variation of @ref step_15 "step-15" that addresses a very
1714 different part of the problem: Line search and whether it is
1715 necessary to re-build the Newton matrix in every nonlinear
1716 iteration. Given that the results above show that using automatic
1717 differentiation comes at a cost, the techniques in @ref step_77 "step-77" have the
1718 potential to offset these costs to some degree. It would therefore
1719 be quite interesting to combine the current program with the
1720 techniques in @ref step_77 "step-77". For production codes, this would certainly be
1721 the way to go.
1722 *
1723 *
1724<a name="step_72-PlainProg"></a>
1725<h1> The plain program</h1>
1726@include "step-72.cc"
1727*/
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
const unsigned int dofs_per_cell
Definition fe_data.h:436
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)
static void initialize(const std::string &filename="", const std::string &output_filename="", const ParameterHandler::OutputStyle output_style_for_output_filename=ParameterHandler::Short, ParameterHandler &prm=ParameterAcceptor::prm, const ParameterHandler::OutputStyle output_style_for_filename=ParameterHandler::DefaultStyle)
void initialize(const MatrixType &A, const AdditionalData &parameters=AdditionalData())
@ wall_times
Definition timer.h:651
#define AssertIndexRange(index, range)
#define AssertThrow(cond, exc)
void mesh_loop(const CellIteratorType &begin, const CellIteratorType &end, const CellWorkerFunctionType &cell_worker, const CopierType &copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const AssembleFlags flags=assemble_own_cells, const BoundaryWorkerFunctionType &boundary_worker=BoundaryWorkerFunctionType(), const FaceWorkerFunctionType &face_worker=FaceWorkerFunctionType(), const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
Definition mesh_loop.h:281
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_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt K
void hyper_ball(Triangulation< dim, spacedim > &tria, const Point< spacedim > &center={}, const double radius=1., const bool attach_spherical_manifold_on_boundary_cells=false)
void refine_and_coarsen_fixed_number(Triangulation< dim, spacedim > &triangulation, const Vector< Number > &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
constexpr char U
constexpr types::blas_int zero
constexpr types::blas_int one
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 cell_residual(Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< Tensor< 1, dim > > &input, const ArrayView< const std::vector< double > > &velocity, double factor=1.)
Definition advection.h:130
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:193
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
Tensor< 2, dim, Number > F(const Tensor< 2, dim, Number > &Grad_u)
constexpr ReturnType< rank, T >::value_type & extract(T &t, const ArrayType &indices)
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)
void copy(const T *begin, const T *end, U *dest)
int(& functions)(const void *v1, const void *v2)
constexpr double PI
Definition numbers.h:239
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
Definition types.h:32
unsigned int boundary_id
Definition types.h:161