deal.II version GIT relicensing-2171-g98af52917d 2024-11-23 14: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
parallel_in_time.h
Go to the documentation of this file.
1
316 *  
317 *   #include "BraidFuncs.hh"
318 *  
319 * @endcode
320 *
321 * This advances the solution forward by one time step.
322 * First some data is collected from the status struct,
323 * namely the start and stop time and the current timestep
324 * number. The timestep size @f$\Delta t@f$ is calculated,
325 * and the step function from the HeatEquation is used to
326 * advance the solution.
327 *
328 * @code
329 *   int my_Step(braid_App app,
330 *   braid_Vector ustop,
331 *   braid_Vector fstop,
332 *   braid_Vector u,
333 *   braid_StepStatus status)
334 *   {
335 *   UNUSED(ustop);
336 *   UNUSED(fstop);
337 *   double tstart; /* current time */
338 *   double tstop; /* evolve to this time*/
339 *   int level;
340 *   double deltaT;
341 *  
342 *   int index;
343 *   braid_StepStatusGetLevel(status, &level);
344 *   braid_StepStatusGetTstartTstop(status, &tstart, &tstop);
345 *   braid_StepStatusGetTIndex(status, &index);
346 *  
347 *   deltaT = tstop - tstart;
348 *  
349 *   ::Vector<double>& solution = u->data;
350 *  
351 *   HeatEquation<2>& heateq = app->eq;
352 *  
353 *   heateq.step(solution, deltaT, tstart, index);
354 *  
355 *   return 0;
356 *   }
357 *  
358 *  
359 * @endcode
360 *
361 * In this function we initialize a vector at an arbitrary time.
362 * At this point we don't know anything about what the solution
363 * looks like, and we can really initialize to anything, so in
364 * this case use reinit to initialize the memory and set the
365 * values to zero.
366 *
367 * @code
368 *   int
369 *   my_Init(braid_App app,
370 *   double t,
371 *   braid_Vector *u_ptr)
372 *   {
373 *   my_Vector *u = new(my_Vector);
374 *   int size = app->eq.size();
375 *   u->data.reinit(size);
376 *  
377 *   app->eq.initialize(t, u->data);
378 *  
379 *   *u_ptr = u;
380 *  
381 *   return 0;
382 *   }
383 *  
384 * @endcode
385 *
386 * Here we need to copy the vector u into the vector v. We do this
387 * by allocating a new vector, then reinitializing the deal.ii
388 * vector to the correct size. The deal.ii reinitialization sets
389 * every value to zero, so next we need to iterate over the vector
390 * u and copy the values to the new vector v.
391 *
392 * @code
393 *   int
394 *   my_Clone(braid_App app,
395 *   braid_Vector u,
396 *   braid_Vector *v_ptr)
397 *   {
398 *   UNUSED(app);
399 *   my_Vector *v = new(my_Vector);
400 *   int size = u->data.size();
401 *   v->data.reinit(size);
402 *   for(size_t i=0, end=v->data.size(); i != end; ++i)
403 *   {
404 *   v->data[i] = u->data[i];
405 *   }
406 *   *v_ptr = v;
407 *  
408 *   return 0;
409 *   }
410 *  
411 * @endcode
412 *
413 * Here we need to free the memory used by vector u. This is
414 * pretty simple since the deal.ii vector is stored inside the
415 * XBraid vector, so we just delete the XBraid vector u and it
416 * puts the deal.ii vector out of scope and releases its memory.
417 *
418 * @code
419 *   int
420 *   my_Free(braid_App app,
421 *   braid_Vector u)
422 *   {
423 *   UNUSED(app);
424 *   delete u;
425 *  
426 *   return 0;
427 *   }
428 *  
429 * @endcode
430 *
431 * This is to perform an axpy type operation. That is to say we
432 * do @f$y = \alpha x + \beta y@f$. Fortunately deal.ii already has
433 * this operation built in to its vector class, so we get the
434 * reference to the vector y and call the sadd method.
435 *
436 * @code
437 *   int my_Sum(braid_App app,
438 *   double alpha,
439 *   braid_Vector x,
440 *   double beta,
441 *   braid_Vector y)
442 *   {
443 *   UNUSED(app);
444 *   Vector<double>& vec = y->data;
445 *   vec.sadd(beta, alpha, x->data);
446 *  
447 *   return 0;
448 *   }
449 *  
450 * @endcode
451 *
452 * This calculates the spatial norm using the l2 norm. According
453 * to XBraid, this could be just about any spatial norm but we'll
454 * keep it simple and used deal.ii vector's built in l2_norm method.
455 *
456 * @code
457 *   int
458 *   my_SpatialNorm(braid_App app,
459 *   braid_Vector u,
460 *   double *norm_ptr)
461 *   {
462 *   UNUSED(app);
463 *   double dot = 0.0;
464 *   dot = u->data.l2_norm();
465 *   *norm_ptr = dot;
466 *  
467 *   return 0;
468 *   }
469 *  
470 * @endcode
471 *
472 * This function is called at various points depending on the access
473 * level specified when configuring the XBraid struct. This function
474 * is used to print out data during the run time, such as plots of the
475 * data. The status struct contains a ton of information about the
476 * simulation run. Here we get the current time and timestep number.
477 * The output_results function is called to plot the solution data.
478 * If the method of manufactured solutions is being used, then the
479 * error of this time step is computed and processed.
480 *
481 * @code
482 *   int
483 *   my_Access(braid_App app,
484 *   braid_Vector u,
485 *   braid_AccessStatus astatus)
486 *   {
487 *   double t;
488 *   int index;
489 *  
490 *   braid_AccessStatusGetT(astatus, &t);
491 *   braid_AccessStatusGetTIndex(astatus, &index);
492 *  
493 *   app->eq.output_results(index, t, u->data);
494 *  
495 *   #if DO_MFG
496 *   if(index == app->final_step)
497 *   {
498 *   app->eq.process_solution(t, index, u->data);
499 *   }
500 *   #endif
501 *  
502 *   return 0;
503 *   }
504 *  
505 * @endcode
506 *
507 * This calculates the size of buffer needed to pack the solution
508 * data into a linear buffer for transfer to another processor via
509 * MPI. We query the size of the data from the HeatEquation class
510 * and return the buffer size.
511 *
512 * @code
513 *   int
514 *   my_BufSize(braid_App app,
515 *   int *size_ptr,
516 *   braid_BufferStatus bstatus)
517 *   {
518 *   UNUSED(bstatus);
519 *   int size = app->eq.size();
520 *   *size_ptr = (size+1)*sizeof(double);
521 *  
522 *   return 0;
523 *   }
524 *  
525 * @endcode
526 *
527 * This function packs a linear buffer with data so that the buffer
528 * may be sent to another processor via MPI. The buffer is cast to
529 * a type we can work with. The first element of the buffer is the
530 * size of the buffer. Then we iterate over solution vector u and
531 * fill the buffer with our solution data. Finally we tell XBraid
532 * how much data we wrote.
533 *
534 * @code
535 *   int
536 *   my_BufPack(braid_App app,
537 *   braid_Vector u,
538 *   void *buffer,
539 *   braid_BufferStatus bstatus)
540 *   {
541 *  
542 *   UNUSED(app);
543 *   double *dbuffer = (double*)buffer;
544 *   int size = u->data.size();
545 *   dbuffer[0] = size;
546 *   for(int i=0; i != size; ++i)
547 *   {
548 *   dbuffer[i+1] = (u->data)[i];
549 *   }
550 *   braid_BufferStatusSetSize(bstatus, (size+1)*sizeof(double));
551 *  
552 *   return 0;
553 *   }
554 *  
555 * @endcode
556 *
557 * This function unpacks a buffer that was received from a different
558 * processor via MPI. The size of the buffer is read from the first
559 * element, then we iterate over the size of the buffer and fill
560 * the values of solution vector u with the data in the buffer.
561 *
562 * @code
563 *   int
564 *   my_BufUnpack(braid_App app,
565 *   void *buffer,
566 *   braid_Vector *u_ptr,
567 *   braid_BufferStatus bstatus)
568 *   {
569 *   UNUSED(app);
570 *   UNUSED(bstatus);
571 *  
572 *   my_Vector *u = NULL;
573 *   double *dbuffer = (double*)buffer;
574 *   int size = static_cast<int>(dbuffer[0]);
575 *   u = new(my_Vector);
576 *   u->data.reinit(size);
577 *  
578 *   for(int i = 0; i != size; ++i)
579 *   {
580 *   (u->data)[i] = dbuffer[i+1];
581 *   }
582 *   *u_ptr = u;
583 *  
584 *   return 0;
585 *   }
586 * @endcode
587
588
589<a name="ann-src/BraidFuncs.hh"></a>
590<h1>Annotated version of src/BraidFuncs.hh</h1>
591 *
592 *
593 *
594 *
595 * @code
596 *   /* -----------------------------------------------------------------------------
597 *   *
598 *   * SPDX-License-Identifier: LGPL-2.1-or-later
599 *   * Copyright (C) 2018 by Joshua Christopher
600 *   *
601 *   * This file is part of the deal.II code gallery.
602 *   *
603 *   * -----------------------------------------------------------------------------
604 *   */
605 *  
606 *   #ifndef _BRAIDFUNCS_H_
607 *   #define _BRAIDFUNCS_H_
608 *  
609 *   /**
610 *   * \file BraidFuncs.cc
611 *   * \brief Contains the implementation of the mandatory X-Braid functions
612 *   *
613 *   * X-Braid mandates several functions in order to drive the solution.
614 *   * This file contains the implementation of said mandatory functions.
615 *   * See the X-Braid documentation for more information.
616 *   * There are several functions that are optional in X-Braid that may
617 *   * or may not be implemented in here.
618 *   *
619 *   */
620 *  
621 *  
622 *   /*-------- Third Party --------*/
623 *   #include <deal.II/numerics/vector_tools.h>
624 *  
625 *   #include <braid.h>
626 *   #include <braid_test.h>
627 *  
628 *   /*-------- Project --------*/
629 *   #include "HeatEquation.hh"
630 *  
631 * @endcode
632 *
633 * This struct contains all data that changes with time. For now
634 * this is just the solution data. When doing AMR this should
635 * probably include the triangulization, the sparsity pattern,
636 * constraints, etc.
637 *
638 * @code
639 *   /**
640 *   * \brief Struct that contains the deal.ii vector.
641 *   */
642 *   typedef struct _braid_Vector_struct
643 *   {
644 *   ::Vector<double> data;
645 *   } my_Vector;
646 *  
647 * @endcode
648 *
649 * This struct contains all the data that is unchanging with time.
650 *
651 * @code
652 *   /**
653 *   * \brief Struct that contains the HeatEquation and final
654 *   * time step number.
655 *   */
656 *   typedef struct _braid_App_struct
657 *   {
658 *   HeatEquation<2> eq;
659 *   int final_step;
660 *   } my_App;
661 *  
662 *  
663 *   /**
664 *   * @brief my_Step - Takes a step in time, advancing the u vector
665 *   *
666 *   * @param app - The braid app struct
667 *   * @param ustop - The solution data at the end of this time step
668 *   * @param fstop - RHS data (such as forcing function?)
669 *   * @param u - The solution data at the beginning of this time step
670 *   * @param status - Status structure that contains various info of this time
671 *   *
672 *   * @return Success (0) or failure (1)
673 *   **/
674 *   int my_Step(braid_App app,
675 *   braid_Vector ustop,
676 *   braid_Vector fstop,
677 *   braid_Vector u,
678 *   braid_StepStatus status);
679 *  
680 *  
681 *   /**
682 *   * @brief my_Init - Initializes a solution data at the given time
683 *   * For now, initializes the solution to zero no matter what time we are at
684 *   *
685 *   * @param app - The braid app struct containing user data
686 *   * @param t - Time at which the solution is initialized
687 *   * @param u_ptr - The solution data that needs to be filled
688 *   *
689 *   * @return Success (0) or failure (1)
690 *   **/
691 *   int
692 *   my_Init(braid_App app,
693 *   double t,
694 *   braid_Vector *u_ptr);
695 *  
696 *  
697 *   /**
698 *   * @brief my_Clone - Clones a vector into a new vector
699 *   *
700 *   * @param app - The braid app struct containing user data
701 *   * @param u - The existing vector containing data
702 *   * @param v_ptr - The empty vector that needs to be filled
703 *   *
704 *   * @return Success (0) or failure (1)
705 *   **/
706 *   int
707 *   my_Clone(braid_App app,
708 *   braid_Vector u,
709 *   braid_Vector *v_ptr);
710 *  
711 *  
712 *   /**
713 *   * @brief my_Free - Deletes a vector
714 *   *
715 *   * @param app - The braid app struct containing user data
716 *   * @param u - The vector that needs to be deleted
717 *   *
718 *   * @return Success (0) or failure (1)
719 *   **/
720 *   int
721 *   my_Free(braid_App app,
722 *   braid_Vector u);
723 *  
724 *  
725 *   /**
726 *   * @brief my_Sum - Sums two vectors in an AXPY operation
727 *   * The operation is y = alpha*x + beta*y
728 *   *
729 *   * @param app - The braid app struct containing user data
730 *   * @param alpha - The coefficient in front of x
731 *   * @param x - A vector that is multiplied by alpha then added to y
732 *   * @param beta - The coefficient of y
733 *   * @param y - A vector that is multiplied by beta then summed with x
734 *   *
735 *   * @return Success (0) or failure (1)
736 *   **/
737 *   int
738 *   my_Sum(braid_App app,
739 *   double alpha,
740 *   braid_Vector x,
741 *   double beta,
742 *   braid_Vector y);
743 *  
744 *   /**
745 *   * \brief Returns the spatial norm of the provided vector
746 *   *
747 *   * Calculates and returns the spatial norm of the provided vector.
748 *   * Interestingly enough, X-Braid does not specify a particular norm.
749 *   * to keep things simple, we implement the Euclidean norm.
750 *   *
751 *   * \param app - The braid app struct containing user data
752 *   * \param u - The vector we need to take the norm of
753 *   * \param norm_ptr - Pointer to the norm that was calculated, need to modify this
754 *   * \return Success (0) or failure (1)
755 *   */
756 *   int
757 *   my_SpatialNorm(braid_App app,
758 *   braid_Vector u,
759 *   double *norm_ptr);
760 *  
761 *   /**
762 *   * \brief Allows the user to output details
763 *   *
764 *   * The Access function is called at various points to allow the user to output
765 *   * information to the screen or to files.
766 *   * The astatus parameter provides various information about the simulation,
767 *   * see the XBraid documentation for details on what information you can get.
768 *   * Example information is what the current timestep number and current time is.
769 *   * If the access level (in parallel_in_time.cc) is set to 0, this function is
770 *   * never called.
771 *   * If the access level is set to 1, the function is called after the last
772 *   * XBraid cycle.
773 *   * If the access level is set to 2, it is called every XBraid cycle.
774 *   *
775 *   * \param app - The braid app struct containing user data
776 *   * \param u - The vector containing the data at the status provided
777 *   * \param astatus - The Braid status structure
778 *   * \return Success (0) or failure (1)
779 *   */
780 *   int
781 *   my_Access(braid_App app,
782 *   braid_Vector u,
783 *   braid_AccessStatus astatus);
784 *  
785 *   /**
786 *   * \brief Calculates the size of a buffer for MPI data transfer
787 *   *
788 *   * Calculates the size of the buffer that is needed to transfer
789 *   * a solution vector to another processor.
790 *   * The bstatus parameter provides various information on the
791 *   * simulation, see the XBraid documentation for all possible
792 *   * fields.
793 *   *
794 *   * \param app - The braid app struct containing user data
795 *   * \param size_ptr A pointer to the calculated size
796 *   * \param bstatus The XBraid status structure
797 *   * \return Success (0) or failure (1)
798 *   */
799 *   int
800 *   my_BufSize(braid_App app,
801 *   int *size_ptr,
802 *   braid_BufferStatus bstatus);
803 *  
804 *   /**
805 *   * \brief Linearizes a vector to be sent to another processor
806 *   *
807 *   * Linearizes (packs) a data buffer with the contents of
808 *   * some solution state u.
809 *   *
810 *   * \param app - The braid app struct containing user data
811 *   * \param u The vector that must be packed into buffer
812 *   * \param buffer The buffer that must be filled with u
813 *   * \param bstatus The XBraid status structure
814 *   * \return Success (0) or failure (1)
815 *   */
816 *   int
817 *   my_BufPack(braid_App app,
818 *   braid_Vector u,
819 *   void *buffer,
820 *   braid_BufferStatus bstatus);
821 *  
822 *   /**
823 *   * \brief Unpacks a vector that was sent from another processor
824 *   *
825 *   * Unpacks a linear data buffer into the vector pointed to by
826 *   * u_ptr.
827 *   *
828 *   * \param app - The braid app struct containing user data
829 *   * \param buffer The buffer that must be unpacked
830 *   * \param u_ptr The pointer to the vector that is filled
831 *   * \param bstatus The XBraid status structure
832 *   * \return Success (0) or failure (1)
833 *   */
834 *   int
835 *   my_BufUnpack(braid_App app,
836 *   void *buffer,
837 *   braid_Vector *u_ptr,
838 *   braid_BufferStatus bstatus);
839 *  
840 *   #endif // _BRAIDFUNCS_H_
841 * @endcode
842
843
844<a name="ann-src/HeatEquation.hh"></a>
845<h1>Annotated version of src/HeatEquation.hh</h1>
846 *
847 *
848 *
849 *
850 * @code
851 *   /* -----------------------------------------------------------------------------
852 *   *
853 *   * SPDX-License-Identifier: LGPL-2.1-or-later
854 *   * Copyright (C) 2018 by Joshua Christopher
855 *   *
856 *   * This file is part of the deal.II code gallery.
857 *   *
858 *   * -----------------------------------------------------------------------------
859 *   */
860 *  
861 *   #ifndef _HEATEQUATION_H_
862 *   #define _HEATEQUATION_H_
863 *  
864 *   #include <deal.II/base/utilities.h>
865 *   #include <deal.II/base/quadrature_lib.h>
866 *   #include <deal.II/base/function.h>
867 *   #include <deal.II/base/logstream.h>
868 *   #include <deal.II/lac/vector.h>
869 *   #include <deal.II/lac/full_matrix.h>
870 *   #include <deal.II/lac/dynamic_sparsity_pattern.h>
871 *   #include <deal.II/lac/sparse_matrix.h>
872 *   #include <deal.II/lac/solver_cg.h>
873 *   #include <deal.II/lac/precondition.h>
874 *   #include <deal.II/lac/affine_constraints.h>
875 *   #include <deal.II/grid/tria.h>
876 *   #include <deal.II/grid/grid_generator.h>
877 *   #include <deal.II/grid/grid_refinement.h>
878 *   #include <deal.II/grid/grid_out.h>
879 *   #include <deal.II/grid/tria_accessor.h>
880 *   #include <deal.II/grid/tria_iterator.h>
881 *   #include <deal.II/dofs/dof_handler.h>
882 *   #include <deal.II/dofs/dof_accessor.h>
883 *   #include <deal.II/dofs/dof_tools.h>
884 *   #include <deal.II/fe/fe_q.h>
885 *   #include <deal.II/fe/fe_values.h>
886 *   #include <deal.II/numerics/data_out.h>
887 *   #include <deal.II/numerics/vector_tools.h>
888 *   #include <deal.II/numerics/error_estimator.h>
889 *   #include <deal.II/numerics/solution_transfer.h>
890 *   #include <deal.II/numerics/matrix_tools.h>
891 *   #include <deal.II/base/convergence_table.h>
892 *  
893 *   #include <fstream>
894 *  
895 *   using namespace dealii;
896 *  
897 * @endcode
898 *
899 * The HeatEquation class is describes the finite element
900 * solver for the heat equation. It contains all the functions
901 * needed to define the problem domain and advance the solution
902 * in time.
903 *
904 * @code
905 *   template <int dim>
906 *   class HeatEquation
907 *   {
908 *   public:
909 *   HeatEquation();
910 *   void define();
911 *   void step(Vector<double>& braid_data,
912 *   double deltaT,
913 *   double a_time,
914 *   int a_time_idx);
915 *  
916 *   int size() const; /// Returns the size of the solution vector
917 *  
918 *   void output_results(int a_time_idx,
919 *   double a_time,
920 *   Vector<double>& a_solution) const;
921 *  
922 *   void initialize(double a_time,
923 *   Vector<double>& a_vector) const;
924 *  
925 *   void process_solution(double a_time,
926 *   int a_index,
927 *   const Vector<double>& a_vector);
928 *  
929 *   private:
930 *   void setup_system();
931 *   void solve_time_step(Vector<double>& a_solution);
932 *  
933 *   Triangulation<dim> triangulation;
934 *   FE_Q<dim> fe;
935 *   DoFHandler<dim> dof_handler;
936 *  
937 *   AffineConstraints<double> constraints;
938 *  
939 *   SparsityPattern sparsity_pattern;
940 *   SparseMatrix<double> mass_matrix;
941 *   SparseMatrix<double> laplace_matrix;
942 *   SparseMatrix<double> system_matrix;
943 *  
944 *   Vector<double> system_rhs;
945 *  
946 *   std::ofstream myfile;
947 *  
948 *   const double theta;
949 *  
950 * @endcode
951 *
952 * These were originally in the run() function but because
953 * I am splitting the run() function up into define and step
954 * they need to become member data
955 *
956 * @code
957 *   Vector<double> tmp;
958 *   Vector<double> forcing_terms;
959 *  
960 *   ConvergenceTable convergence_table;
961 *   };
962 *  
963 * @endcode
964 *
965 * The RightHandSide class describes the RHS of the governing
966 * equations. In this case, it is the forcing function.
967 *
968 * @code
969 *   template <int dim>
970 *   class RightHandSide : public Function<dim>
971 *   {
972 *   public:
973 *   RightHandSide ()
974 *   :
975 *   Function<dim>(),
976 *   period (0.2)
977 *   {}
978 *  
979 *   virtual double value (const Point<dim> &p,
980 *   const unsigned int component = 0) const override;
981 *  
982 *   private:
983 *   const double period;
984 *   };
985 *  
986 * @endcode
987 *
988 * The BoundaryValues class describes the boundary conditions
989 * of the governing equations.
990 *
991 * @code
992 *   template <int dim>
993 *   class BoundaryValues : public Function<dim>
994 *   {
995 *   public:
996 *   virtual double value (const Point<dim> &p,
997 *   const unsigned int component = 0) const override;
998 *   };
999 *  
1000 * @endcode
1001 *
1002 * The RightHandSideMFG class describes the right hand side
1003 * function when doing the method of manufactured solutions.
1004 *
1005 * @code
1006 *   template <int dim>
1007 *   class RightHandSideMFG : public Function<dim>
1008 *   {
1009 *   public:
1010 *   virtual double value (const Point<dim> &p,
1011 *   const unsigned int component = 0) const override;
1012 *   };
1013 *  
1014 * @endcode
1015 *
1016 * The InitialValuesMFG class describes the initial values
1017 * when doing the method of manufactured solutions.
1018 *
1019 * @code
1020 *   template <int dim>
1021 *   class InitialValuesMFG : public Function<dim>
1022 *   {
1023 *   public:
1024 *   virtual double value (const Point<dim> &p,
1025 *   const unsigned int component = 0) const override;
1026 *   };
1027 *  
1028 * @endcode
1029 *
1030 * Provides the exact value for the manufactured solution. This
1031 * is used for the boundary conditions as well.
1032 *
1033 * @code
1034 *   template <int dim>
1035 *   class ExactValuesMFG : public Function<dim>
1036 *   {
1037 *   public:
1038 *   /**
1039 *   * \brief Computes the value at the given point and member data time
1040 *   *
1041 *   * Computes the exact value of the manufactured solution at point p and
1042 *   * the member data time. See the class documentation and the design doc
1043 *   * for details on what the exact solution is.
1044 *   *
1045 *   * \param p The point that the exact solution is computed at
1046 *   * \param component The component of the exact solution (always 0 for now)
1047 *   * \return double The exact value that was computed
1048 *   */
1049 *   virtual double value (const Point<dim> &p,
1050 *   const unsigned int component = 0) const override;
1051 *  
1052 *   /**
1053 *   * \brief Computes the gradient of the exact solution at the given point
1054 *   *
1055 *   * Computes the gradient of the exact/manufactured solution value at
1056 *   * point p and member data time. See the design doc for details on
1057 *   * what the gradient of the exact solution is
1058 *   *
1059 *   * \param p The point that the gradient is calculated at
1060 *   * \param component The component of the system of equations this gradient is for
1061 *   * \return Tensor<1,dim> A rank 1 tensor that contains the gradient
1062 *   * in each spatial dimension
1063 *   */
1064 *   virtual Tensor<1,dim> gradient (const Point<dim> &p,
1065 *   const unsigned int component = 0) const override;
1066 *   };
1067 *  
1068 *  
1069 *   #include "HeatEquationImplem.hh"
1070 *  
1071 *   #endif // _HEATEQUATION_H_
1072 * @endcode
1073
1074
1075<a name="ann-src/HeatEquationImplem.hh"></a>
1076<h1>Annotated version of src/HeatEquationImplem.hh</h1>
1077 *
1078 *
1079 *
1080 *
1081 * @code
1082 *   /* -----------------------------------------------------------------------------
1083 *   *
1084 *   * SPDX-License-Identifier: LGPL-2.1-or-later
1085 *   * Copyright (C) 2018 by Joshua Christopher
1086 *   *
1087 *   * This file is part of the deal.II code gallery.
1088 *   *
1089 *   * -----------------------------------------------------------------------------
1090 *   */
1091 *  
1092 *   #include "Utilities.hh"
1093 *  
1094 *   #include <iomanip>
1095 *   #include <math.h>
1096 *  
1097 * @endcode
1098 *
1099 * Calculates the forcing function for the RightHandSide. See the
1100 * documentation for the math.
1101 *
1102 * @code
1103 *   template <int dim>
1104 *   double RightHandSide<dim>::value (const Point<dim> &p,
1105 *   const unsigned int component) const
1106 *   {
1107 *   (void) component;
1108 *   Assert (component == 0, ExcIndexRange(component, 0, 1));
1109 *   Assert (dim == 2, ExcNotImplemented());
1110 *  
1111 *   double time = this->get_time();
1112 *  
1113 *   if ((p[0] > 0.5) && (p[1] > -0.5))
1114 *   {
1115 *   return std::exp(-0.5*(time-0.125)*(time-0.125)/(0.005));
1116 *   }
1117 *   else if ((p[0] > -0.5) && (p[1] > 0.5))
1118 *   {
1119 *   return std::exp(-0.5*(time-0.375)*(time-0.375)/(0.005));
1120 *   }
1121 *   else
1122 *   {
1123 *   return 0;
1124 *   }
1125 *  
1126 *   return 0; // No forcing function
1127 *   }
1128 *  
1129 * @endcode
1130 *
1131 * Calculates the forcing function for the method of manufactured
1132 * solutions. See the documentation for the math.
1133 *
1134 * @code
1135 *   template <int dim>
1136 *   double RightHandSideMFG<dim>::value (const Point<dim> &p,
1137 *   const unsigned int component) const
1138 *   {
1139 *   (void) component;
1140 *   Assert (component == 0, ExcIndexRange(component, 0, 1));
1141 *   Assert (dim == 2, ExcNotImplemented());
1142 *  
1143 *   double time = this->get_time();
1144 *  
1145 *   double pi = numbers::PI;
1146 *   return 4*pi*pi*std::exp(-4*pi*pi*time)*std::cos(2*pi*p[0])*std::cos(2*pi*p[1]);
1147 *   }
1148 *  
1149 * @endcode
1150 *
1151 * Calculates the boundary conditions, essentially zero everywhere.
1152 *
1153 * @code
1154 *   template <int dim>
1155 *   double BoundaryValues<dim>::value (const Point<dim> &p,
1156 *   const unsigned int component) const
1157 *   {
1158 *   UNUSED(p);
1159 *   (void) component;
1160 *   Assert (component == 0, ExcIndexRange(component, 0, 1));
1161 *   return 0;
1162 *   }
1163 *  
1164 * @endcode
1165 *
1166 * Calculates the exact solution (and thus also boundary conditions)
1167 * for the method of manufactured solutions.
1168 *
1169 * @code
1170 *   template <int dim>
1171 *   double ExactValuesMFG<dim>::value (const Point<dim> &p,
1172 *   const unsigned int component) const
1173 *   {
1174 *   (void) component;
1175 *   Assert (component == 0, ExcIndexRange(component, 0, 1));
1176 *  
1177 *   double time = this->get_time();
1178 *   const double pi = numbers::PI;
1179 *  
1180 *   return std::exp(-4*pi*pi*time)*std::cos(2*pi*p[0])*std::cos(2*pi*p[1]);
1181 *   }
1182 *  
1183 * @endcode
1184 *
1185 * Calculates the gradient of the exact solution for the method of manufactured
1186 * solutions. See the documentation for the math.
1187 *
1188 * @code
1189 *   template <int dim>
1190 *   Tensor<1,dim> ExactValuesMFG<dim>::gradient (const Point<dim> &p,
1191 *   const unsigned int) const
1192 *   {
1193 *   Assert (dim == 2, ExcNotImplemented());
1194 *  
1195 *   Tensor<1,dim> return_value;
1196 *   const double pi = numbers::PI;
1197 *   double time = this->get_time();
1198 *   return_value[0] = -2*pi*std::exp(-4*pi*pi*time)*std::cos(2*pi*p[1])*std::sin(2*pi*p[0]);
1199 *   return_value[1] = -2*pi*std::exp(-4*pi*pi*time)*std::cos(2*pi*p[0])*std::sin(2*pi*p[1]);
1200 *   return return_value;
1201 *   }
1202 *  
1203 * @endcode
1204 *
1205 * Calculates the initial values for the method of manufactured solutions.
1206 * See the documentation for the math.
1207 *
1208 * @code
1209 *   template <int dim>
1210 *   double InitialValuesMFG<dim>::value (const Point<dim> &p,
1211 *   const unsigned int component) const
1212 *   {
1213 *   (void) component;
1214 *   Assert (component == 0, ExcIndexRange(component, 0, 1));
1215 *   const double pi = numbers::PI;
1216 *  
1217 *   return std::cos(2*pi*p[0])*std::cos(2*pi*p[1]);
1218 *   }
1219 *  
1220 *   template <int dim>
1221 *   HeatEquation<dim>::HeatEquation ()
1222 *   :
1223 *   fe(1),
1224 *   dof_handler(triangulation),
1225 *   theta(0.5)
1226 *   {
1227 *   }
1228 *  
1229 *   template <int dim>
1230 *   void HeatEquation<dim>::initialize(double a_time,
1231 *   Vector<double>& a_vector) const
1232 *   {
1233 *   #if DO_MFG
1234 * @endcode
1235 *
1236 * We only initialize values in the manufactured solution case
1237 *
1238 * @code
1239 *   InitialValuesMFG<dim> iv_function;
1240 *   iv_function.set_time(a_time);
1241 *   VectorTools::project (dof_handler, constraints,
1242 *   QGauss<dim>(fe.degree+1), iv_function,
1243 *   a_vector);
1244 *   #else
1245 *   UNUSED(a_time);
1246 *   UNUSED(a_vector);
1247 *   #endif // DO_MFG
1248 * @endcode
1249 *
1250 * If not the MFG solution case, a_vector is already zero'd so do nothing
1251 *
1252 * @code
1253 *   }
1254 *  
1255 *   template <int dim>
1256 *   void HeatEquation<dim>::setup_system()
1257 *   {
1258 *   dof_handler.distribute_dofs(fe);
1259 *  
1260 *   constraints.clear ();
1262 *   constraints);
1263 *   constraints.close();
1264 *  
1265 *   DynamicSparsityPattern dsp(dof_handler.n_dofs());
1266 *   DoFTools::make_sparsity_pattern(dof_handler,
1267 *   dsp,
1268 *   constraints,
1269 *   /*keep_constrained_dofs = */ true);
1270 *   sparsity_pattern.copy_from(dsp);
1271 *  
1272 *   mass_matrix.reinit(sparsity_pattern);
1273 *   laplace_matrix.reinit(sparsity_pattern);
1274 *   system_matrix.reinit(sparsity_pattern);
1275 *  
1276 *   MatrixCreator::create_mass_matrix(dof_handler,
1277 *   QGauss<dim>(fe.degree+1),
1278 *   mass_matrix);
1279 *   MatrixCreator::create_laplace_matrix(dof_handler,
1280 *   QGauss<dim>(fe.degree+1),
1281 *   laplace_matrix);
1282 *  
1283 *   system_rhs.reinit(dof_handler.n_dofs());
1284 *   }
1285 *  
1286 *  
1287 *   template <int dim>
1288 *   void HeatEquation<dim>::solve_time_step(Vector<double>& a_solution)
1289 *   {
1290 *   SolverControl solver_control(1000, 1e-8 * system_rhs.l2_norm());
1291 *   SolverCG<> cg(solver_control);
1292 *  
1293 *   PreconditionSSOR<> preconditioner;
1294 *   preconditioner.initialize(system_matrix, 1.0);
1295 *  
1296 *   cg.solve(system_matrix, a_solution, system_rhs,
1297 *   preconditioner);
1298 *  
1299 *   constraints.distribute(a_solution);
1300 *   }
1301 *  
1302 *  
1303 *  
1304 *   template <int dim>
1305 *   void HeatEquation<dim>::output_results(int a_time_idx,
1306 *   double a_time,
1307 *   Vector<double>& a_solution) const
1308 *   {
1309 *  
1310 *   DataOutBase::VtkFlags vtk_flags;
1311 *   vtk_flags.time = a_time;
1312 *   vtk_flags.cycle = a_time_idx;
1313 *  
1314 *   DataOut<dim> data_out;
1315 *   data_out.set_flags(vtk_flags);
1316 *  
1317 *   data_out.attach_dof_handler(dof_handler);
1318 *   data_out.add_data_vector(a_solution, "U");
1319 *  
1320 *   data_out.build_patches();
1321 *  
1322 *   const std::string filename = "solution-"
1323 *   + Utilities::int_to_string(a_time_idx, 3) +
1324 *   ".vtk";
1325 *   std::ofstream output(filename.c_str());
1326 *   data_out.write_vtk(output);
1327 *   }
1328 *  
1329 * @endcode
1330 *
1331 * We define the geometry here, this is called on each processor
1332 * and doesn't change in time. Once doing AMR, this won't need
1333 * to exist anymore.
1334 *
1335 * @code
1336 *   template <int dim>
1337 *   void HeatEquation<dim>::define()
1338 *   {
1339 *   const unsigned int initial_global_refinement = 6;
1340 *  
1342 *   triangulation.refine_global (initial_global_refinement);
1343 *  
1344 *   setup_system();
1345 *  
1346 *   tmp.reinit (dof_handler.n_dofs());
1347 *   forcing_terms.reinit (dof_handler.n_dofs());
1348 *   }
1349 *  
1350 * @endcode
1351 *
1352 * Here we advance the solution forward in time. This is done
1353 * the same way as in the loop in @ref step_26 "step-26"'s run function.
1354 *
1355 * @code
1356 *   template<int dim>
1357 *   void HeatEquation<dim>::step(Vector<double>& braid_data,
1358 *   double deltaT,
1359 *   double a_time,
1360 *   int a_time_idx)
1361 *   {
1362 *   a_time += deltaT;
1363 *   ++a_time_idx;
1364 *  
1365 *   mass_matrix.vmult(system_rhs, braid_data);
1366 *  
1367 *   laplace_matrix.vmult(tmp, braid_data);
1368 *  
1369 *   system_rhs.add(-(1 - theta) * deltaT, tmp);
1370 *  
1371 *   #if DO_MFG
1372 *   RightHandSideMFG<dim> rhs_function;
1373 *   #else
1374 *   RightHandSide<dim> rhs_function;
1375 *   #endif
1376 *   rhs_function.set_time(a_time);
1377 *   VectorTools::create_right_hand_side(dof_handler,
1378 *   QGauss<dim>(fe.degree+1),
1379 *   rhs_function,
1380 *   tmp);
1381 *  
1382 *   forcing_terms = tmp;
1383 *   forcing_terms *= deltaT * theta;
1384 *  
1385 *   rhs_function.set_time(a_time - deltaT);
1386 *   VectorTools::create_right_hand_side(dof_handler,
1387 *   QGauss<dim>(fe.degree+1),
1388 *   rhs_function,
1389 *   tmp);
1390 *  
1391 *   forcing_terms.add(deltaT * (1 - theta), tmp);
1392 *   system_rhs += forcing_terms;
1393 *  
1394 *   system_matrix.copy_from(mass_matrix);
1395 *   system_matrix.add(theta * deltaT, laplace_matrix);
1396 *  
1397 *   constraints.condense (system_matrix, system_rhs);
1398 *  
1399 *   {
1400 *   #if DO_MFG
1401 * @endcode
1402 *
1403 * If we are doing the method of manufactured solutions
1404 * then we set the boundary conditions to the exact solution.
1405 * Otherwise the boundary conditions are zero.
1406 *
1407 * @code
1408 *   ExactValuesMFG<dim> boundary_values_function;
1409 *   #else
1410 *   BoundaryValues<dim> boundary_values_function;
1411 *   #endif
1412 *   boundary_values_function.set_time(a_time);
1413 *  
1414 *   std::map<types::global_dof_index, double> boundary_values;
1415 *   VectorTools::interpolate_boundary_values(dof_handler,
1416 *   0,
1417 *   boundary_values_function,
1418 *   boundary_values);
1419 *  
1420 *   MatrixTools::apply_boundary_values(boundary_values,
1421 *   system_matrix,
1422 *   braid_data,
1423 *   system_rhs);
1424 *   }
1425 *  
1426 *   solve_time_step(braid_data);
1427 *   }
1428 *  
1429 *   template<int dim>
1430 *   int HeatEquation<dim>::size() const
1431 *   {
1432 *   return dof_handler.n_dofs();
1433 *   }
1434 *  
1435 * @endcode
1436 *
1437 * This function computes the error for the time step when doing
1438 * the method of manufactured solutions. First the exact values
1439 * is calculated, then the difference per cell is computed for
1440 * the various norms, and the error is computed. This is written
1441 * out to a pretty table.
1442 *
1443 * @code
1444 *   template<int dim> void
1445 *   HeatEquation<dim>::process_solution(double a_time,
1446 *   int a_index,
1447 *   const Vector<double>& a_vector)
1448 *   {
1449 * @endcode
1450 *
1451 * Compute the exact value for the manufactured solution case
1452 *
1453 * @code
1454 *   ExactValuesMFG<dim> exact_function;
1455 *   exact_function.set_time(a_time);
1456 *  
1457 *   Vector<double> difference_per_cell (triangulation.n_active_cells());
1458 *   VectorTools::integrate_difference(dof_handler,
1459 *   a_vector,
1460 *   exact_function,
1461 *   difference_per_cell,
1462 *   QGauss<dim>(fe.degree+1),
1463 *   VectorTools::L2_norm);
1464 *  
1465 *   const double L2_error = VectorTools::compute_global_error(triangulation,
1466 *   difference_per_cell,
1467 *   VectorTools::L2_norm);
1468 *  
1469 *   VectorTools::integrate_difference(dof_handler,
1470 *   a_vector,
1471 *   exact_function,
1472 *   difference_per_cell,
1473 *   QGauss<dim>(fe.degree+1),
1474 *   VectorTools::H1_seminorm);
1475 *  
1476 *   const double H1_error = VectorTools::compute_global_error(triangulation,
1477 *   difference_per_cell,
1478 *   VectorTools::H1_seminorm);
1479 *  
1480 *   const QTrapezoid<1> q_trapez;
1481 *   const QIterated<dim> q_iterated (q_trapez, 5);
1482 *   VectorTools::integrate_difference (dof_handler,
1483 *   a_vector,
1484 *   exact_function,
1485 *   difference_per_cell,
1486 *   q_iterated,
1487 *   VectorTools::Linfty_norm);
1488 *   const double Linfty_error = VectorTools::compute_global_error(triangulation,
1489 *   difference_per_cell,
1490 *   VectorTools::Linfty_norm);
1491 *  
1492 *   const unsigned int n_active_cells = triangulation.n_active_cells();
1493 *   const unsigned int n_dofs = dof_handler.n_dofs();
1494 *  
1495 *   pout() << "Cycle " << a_index << ':'
1496 *   << std::endl
1497 *   << " Number of active cells: "
1498 *   << n_active_cells
1499 *   << std::endl
1500 *   << " Number of degrees of freedom: "
1501 *   << n_dofs
1502 *   << std::endl;
1503 *  
1504 *   convergence_table.add_value("cycle", a_index);
1505 *   convergence_table.add_value("cells", n_active_cells);
1506 *   convergence_table.add_value("dofs", n_dofs);
1507 *   convergence_table.add_value("L2", L2_error);
1508 *   convergence_table.add_value("H1", H1_error);
1509 *   convergence_table.add_value("Linfty", Linfty_error);
1510 *  
1511 *   convergence_table.set_precision("L2", 3);
1512 *   convergence_table.set_precision("H1", 3);
1513 *   convergence_table.set_precision("Linfty", 3);
1514 *  
1515 *   convergence_table.set_scientific("L2", true);
1516 *   convergence_table.set_scientific("H1", true);
1517 *   convergence_table.set_scientific("Linfty", true);
1518 *  
1519 *   convergence_table.set_tex_caption("cells", "\\# cells");
1520 *   convergence_table.set_tex_caption("dofs", "\\# dofs");
1521 *   convergence_table.set_tex_caption("L2", "@fL^2@f-error");
1522 *   convergence_table.set_tex_caption("H1", "@fH^1@f-error");
1523 *   convergence_table.set_tex_caption("Linfty", "@fL^\\infty@f-error");
1524 *  
1525 *   convergence_table.set_tex_format("cells", "r");
1526 *   convergence_table.set_tex_format("dofs", "r");
1527 *  
1528 *   std::cout << std::endl;
1529 *   convergence_table.write_text(std::cout);
1530 *  
1531 *   std::ofstream error_table_file("tex-conv-table.tex");
1532 *   convergence_table.write_tex(error_table_file);
1533 *   }
1534 * @endcode
1535
1536
1537<a name="ann-src/Utilities.cc"></a>
1538<h1>Annotated version of src/Utilities.cc</h1>
1539 *
1540 *
1541 *
1542 *
1543 * @code
1544 *   /* -----------------------------------------------------------------------------
1545 *   *
1546 *   * SPDX-License-Identifier: LGPL-2.1-or-later
1547 *   * Copyright (C) 2018 by Joshua Christopher
1548 *   *
1549 *   * This file is part of the deal.II code gallery.
1550 *   *
1551 *   * -----------------------------------------------------------------------------
1552 *   */
1553 *  
1554 *   #include "Utilities.hh"
1555 *  
1556 *   #include <string>
1557 *   #include <fstream>
1558 *  
1559 *   #include <mpi.h>
1560 *  
1561 *   int procID = 0;
1562 *  
1563 *   #ifdef USE_MPI
1564 * @endcode
1565 *
1566 * The shared variables
1567 *
1568
1569 *
1570 *
1571 * @code
1572 *   static std::string s_pout_filename ;
1573 *   static std::string s_pout_basename ;
1574 *   static std::ofstream s_pout ;
1575 *  
1576 *   static bool s_pout_init = false ;
1577 *   static bool s_pout_open = false ;
1578 *  
1579 * @endcode
1580 *
1581 * in parallel, compute the filename give the basename
1582 * [NOTE: dont call this before MPI is initialized.]
1583 *
1584 * @code
1585 *   static void setFileName()
1586 *   {
1587 *   static const size_t ProcnumSize = 1 + 10 + 1 ; //'.' + 10digits + '\0'
1588 *   char procnum[ProcnumSize] ;
1589 *   snprintf( procnum ,ProcnumSize ,".%d" ,procID);
1590 *   s_pout_filename = s_pout_basename + procnum ;
1591 *   }
1592 *  
1593 * @endcode
1594 *
1595 * in parallel, close the file if nec., open it and check for success
1596 *
1597 * @code
1598 *   static void openFile()
1599 *   {
1600 *   if ( s_pout_open )
1601 *   {
1602 *   s_pout.close();
1603 *   }
1604 *   s_pout.open( s_pout_filename.c_str() );
1605 * @endcode
1606 *
1607 * if open() fails, we have problems, but it's better
1608 * to try again later than to make believe it succeeded
1609 *
1610 * @code
1611 *   s_pout_open = (bool)s_pout ;
1612 *   }
1613 *   #endif
1614 *  
1615 *   std::ostream& pout()
1616 *   {
1617 *   #ifdef USE_MPI
1618 * @endcode
1619 *
1620 * the common case is _open == true, which just returns s_pout
1621 *
1622 * @code
1623 *   if ( ! s_pout_open )
1624 *   {
1625 * @endcode
1626 *
1627 * the uncommon cae: the file isn't opened, MPI may not be
1628 * initialized, and the basename may not have been set
1629 *
1630 * @code
1631 *   int flag_i, flag_f;
1632 *   MPI_Initialized(&flag_i);
1633 *   MPI_Finalized(&flag_f);
1634 * @endcode
1635 *
1636 * app hasn't set a basename yet, so set the default
1637 *
1638 * @code
1639 *   if ( ! s_pout_init )
1640 *   {
1641 *   s_pout_basename = "pout" ;
1642 *   s_pout_init = true ;
1643 *   }
1644 * @endcode
1645 *
1646 * if MPI not initialized, we can't open the file so return cout
1647 *
1648 * @code
1649 *   if ( ! flag_i || flag_f)
1650 *   {
1651 *   return std::cout; // MPI hasn't been started yet, or has ended....
1652 *   }
1653 * @endcode
1654 *
1655 * MPI is initialized, so file must not be, so open it
1656 *
1657 * @code
1658 *   setFileName() ;
1659 *   openFile() ;
1660 * @endcode
1661 *
1662 * finally, in case the open failed, return cout
1663 *
1664 * @code
1665 *   if ( ! s_pout_open )
1666 *   {
1667 *   return std::cout ;
1668 *   }
1669 *   }
1670 *   return s_pout ;
1671 *   #else
1672 *   return std::cout;
1673 *   #endif
1674 *   }
1675 * @endcode
1676
1677
1678<a name="ann-src/Utilities.hh"></a>
1679<h1>Annotated version of src/Utilities.hh</h1>
1680 *
1681 *
1682 *
1683 *
1684 * @code
1685 *   /* -----------------------------------------------------------------------------
1686 *   *
1687 *   * SPDX-License-Identifier: LGPL-2.1-or-later
1688 *   * Copyright (C) 2018 by Joshua Christopher
1689 *   *
1690 *   * This file is part of the deal.II code gallery.
1691 *   *
1692 *   * -----------------------------------------------------------------------------
1693 *   */
1694 *  
1695 *   #ifndef _UTILITIES_H_
1696 *   #define _UTILITIES_H_
1697 *  
1698 *   #include <iostream>
1699 *  
1700 * @endcode
1701 *
1702 * This preprocessor macro is used on function arguments
1703 * that are not used in the function. It is used to
1704 * suppress compiler warnings.
1705 *
1706 * @code
1707 *   #define UNUSED(x) (void)(x)
1708 *  
1709 * @endcode
1710 *
1711 * Contains the current MPI processor ID.
1712 *
1713 * @code
1714 *   extern int procID;
1715 *  
1716 * @endcode
1717 *
1718 * Function to return the ostream to write out to. In MPI
1719 * mode it returns a stream to a file named pout.<#> where
1720 * <#> is the procID. This allows the user to write output
1721 * from each processor to a separate file. In serial mode
1722 * (no MPI), it returns the standard output.
1723 *
1724 * @code
1725 *   std::ostream& pout();
1726 *   #endif
1727 * @endcode
1728
1729
1730<a name="ann-src/parallel_in_time.cc"></a>
1731<h1>Annotated version of src/parallel_in_time.cc</h1>
1732 *
1733 *
1734 *
1735 *
1736 * @code
1737 *   /* -----------------------------------------------------------------------------
1738 *   *
1739 *   * SPDX-License-Identifier: LGPL-2.1-or-later
1740 *   * Copyright (C) 2018 by Joshua Christopher
1741 *   *
1742 *   * This file is part of the deal.II code gallery.
1743 *   *
1744 *   * -----------------------------------------------------------------------------
1745 *   *
1746 *   * Author: Joshua Christopher, Colorado State University, 2018
1747 *   */
1748 *  
1749 *   #include "BraidFuncs.hh"
1750 *   #include "HeatEquation.hh"
1751 *   #include "Utilities.hh"
1752 *  
1753 *   #include <fstream>
1754 *   #include <iostream>
1755 *  
1756 *   int main(int argc, char *argv[])
1757 *   {
1758 *   try
1759 *   {
1760 *   using namespace dealii;
1761 *  
1762 *   /* Initialize MPI */
1763 *   MPI_Comm comm; //, comm_x, comm_t;
1764 *   int rank;
1765 *   MPI_Init(&argc, &argv);
1766 *   comm = MPI_COMM_WORLD;
1767 *   MPI_Comm_rank(comm, &rank);
1768 *   procID = rank;
1769 *  
1770 * @endcode
1771 *
1772 * Set up X-Braid
1773 *
1774 * @code
1775 *   /* Initialize Braid */
1776 *   braid_Core core;
1777 *   double tstart = 0.0;
1778 *   double tstop = 0.002;
1779 *   int ntime = 10;
1780 *   my_App *app = new(my_App);
1781 *  
1782 *   braid_Init(MPI_COMM_WORLD, comm, tstart, tstop, ntime, app,
1783 *   my_Step, my_Init, my_Clone, my_Free, my_Sum, my_SpatialNorm,
1784 *   my_Access, my_BufSize, my_BufPack, my_BufUnpack, &core);
1785 *  
1786 *   /* Define XBraid parameters
1787 *   * See -help message forf descriptions */
1788 *   int max_levels = 3;
1789 * @endcode
1790 *
1791 * int nrelax = 1;
1792 * int skip = 0;
1793 *
1794 * @code
1795 *   double tol = 1.e-7;
1796 * @endcode
1797 *
1798 * int cfactor = 2;
1799 *
1800 * @code
1801 *   int max_iter = 5;
1802 * @endcode
1803 *
1804 * int min_coarse = 10;
1805 * int fmg = 0;
1806 * int scoarsen = 0;
1807 * int res = 0;
1808 * int wrapper_tests = 0;
1809 *
1810 * @code
1811 *   int print_level = 1;
1812 *   int access_level = 1;
1813 *   int use_sequential= 0;
1814 *  
1815 *   braid_SetPrintLevel( core, print_level);
1816 *   braid_SetAccessLevel( core, access_level);
1817 *   braid_SetMaxLevels(core, max_levels);
1818 * @endcode
1819 *
1820 * braid_SetMinCoarse( core, min_coarse );
1821 * braid_SetSkip(core, skip);
1822 * braid_SetNRelax(core, -1, nrelax);
1823 *
1824 * @code
1825 *   braid_SetAbsTol(core, tol);
1826 * @endcode
1827 *
1828 * braid_SetCFactor(core, -1, cfactor);
1829 *
1830 * @code
1831 *   braid_SetMaxIter(core, max_iter);
1832 *   braid_SetSeqSoln(core, use_sequential);
1833 *  
1834 *   app->eq.define();
1835 *   app->final_step = ntime;
1836 *  
1837 *   braid_Drive(core);
1838 *  
1839 * @endcode
1840 *
1841 * Free the memory now that we are done
1842 *
1843 * @code
1844 *   braid_Destroy(core);
1845 *  
1846 *   delete app;
1847 *  
1848 * @endcode
1849 *
1850 * Clean up MPI
1851 * MPI_Comm_free(&comm);
1852 *
1853 * @code
1854 *   MPI_Finalize();
1855 *   }
1856 *   catch (std::exception &exc)
1857 *   {
1858 *   std::cerr << std::endl << std::endl
1859 *   << "----------------------------------------------------"
1860 *   << std::endl;
1861 *   std::cerr << "Exception on processing: " << std::endl << exc.what()
1862 *   << std::endl << "Aborting!" << std::endl
1863 *   << "----------------------------------------------------"
1864 *   << std::endl;
1865 *  
1866 *   return 1;
1867 *   }
1868 *   catch (...)
1869 *   {
1870 *   std::cerr << std::endl << std::endl
1871 *   << "----------------------------------------------------"
1872 *   << std::endl;
1873 *   std::cerr << "Unknown exception!" << std::endl << "Aborting!"
1874 *   << std::endl
1875 *   << "----------------------------------------------------"
1876 *   << std::endl;
1877 *   return 1;
1878 *   }
1879 *  
1880 *   return 0;
1881 *   }
1882 *  
1883 * @endcode
1884
1885
1886<a name="ann-test/test_braid.cc"></a>
1887<h1>Annotated version of test/test_braid.cc</h1>
1888 *
1889 *
1890 *
1891 *
1892 * @code
1893 *   /* -----------------------------------------------------------------------------
1894 *   *
1895 *   * SPDX-License-Identifier: LGPL-2.1-or-later
1896 *   * Copyright (C) 2018 by Joshua Christopher
1897 *   *
1898 *   * This file is part of the deal.II code gallery.
1899 *   *
1900 *   * -----------------------------------------------------------------------------
1901 *   */
1902 *  
1903 *   #include "BraidFuncs.hh"
1904 *  
1905 *   #include <braid.h>
1906 *   #include <braid_test.h>
1907 *  
1908 *   #include <iostream>
1909 *  
1910 *   int main(int argc, char** argv)
1911 *   {
1912 *   MPI_Comm comm;
1913 *   int rank;
1914 *   MPI_Init(&argc, &argv);
1915 *   comm = MPI_COMM_WORLD;
1916 *   MPI_Comm_rank(comm, &rank);
1917 *  
1918 *   my_App *app = new(my_App);
1919 *   app->eq.define();
1920 *  
1921 *   double time = 0.2;
1922 *  
1923 *   braid_Int init_access_result = braid_TestInitAccess(app,
1924 *   comm,
1925 *   stdout,
1926 *   time,
1927 *   my_Init,
1928 *   my_Access,
1929 *   my_Free);
1930 *   (void)init_access_result;
1931 *  
1932 *   braid_Int clone_result = braid_TestClone(app,
1933 *   comm,
1934 *   stdout,
1935 *   time,
1936 *   my_Init,
1937 *   my_Access,
1938 *   my_Free,
1939 *   my_Clone);
1940 *   (void)clone_result;
1941 *  
1942 *   braid_Int sum_result = braid_TestSum(app,
1943 *   comm,
1944 *   stdout,
1945 *   time,
1946 *   my_Init,
1947 *   my_Access,
1948 *   my_Free,
1949 *   my_Clone,
1950 *   my_Sum);
1951 *   (void)sum_result;
1952 *  
1953 *   braid_Int norm_result = braid_TestSpatialNorm(app,
1954 *   comm,
1955 *   stdout,
1956 *   time,
1957 *   my_Init,
1958 *   my_Free,
1959 *   my_Clone,
1960 *   my_Sum,
1961 *   my_SpatialNorm);
1962 *   (void)norm_result;
1963 *  
1964 *   braid_Int buf_result = braid_TestBuf(app,
1965 *   comm,
1966 *   stdout,
1967 *   time,
1968 *   my_Init,
1969 *   my_Free,
1970 *   my_Sum,
1971 *   my_SpatialNorm,
1972 *   my_BufSize,
1973 *   my_BufPack,
1974 *   my_BufUnpack);
1975 *   (void)buf_result;
1976 *  
1977 * @endcode
1978 *
1979 * /* Create spatial communicator for wrapper-tests */
1980 * braid_SplitCommworld(&comm, 1, &comm_x, &comm_t);
1981 *
1982 *
1983 * braid_TestAll(app, comm_x, stdout, 0.0, (tstop-tstart)/ntime,
1984 * 2*(tstop-tstart)/ntime, my_Init, my_Free, my_Clone,
1985 * my_Sum, my_SpatialNorm, my_BufSize, my_BufPack,
1986 * my_BufUnpack, my_Coarsen, my_Interp, my_Residual, my_Step);
1987 *
1988
1989 *
1990 *
1991 * @code
1992 *   /* Finalize MPI */
1994 *   }
1995 * @endcode
1996
1997
1998*/
void set_flags(const FlagType &flags)
void initialize(const MatrixType &A, const AdditionalData &parameters=AdditionalData())
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_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)
std::vector< index_type > data
Definition mpi.cc:735
std::size_t size
Definition mpi.cc:734
void hyper_L(Triangulation< dim > &tria, const double left=-1., const double right=1., const bool colorize=false)
void mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
Definition l2.h:57
void create_mass_matrix(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, MatrixType &matrix, const Function< spacedim, typename MatrixType::value_type > *const a=nullptr, const AffineConstraints< typename MatrixType::value_type > &constraints=AffineConstraints< typename MatrixType::value_type >())
void create_laplace_matrix(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, MatrixType &matrix, const Function< spacedim, typename MatrixType::value_type > *const a=nullptr, const AffineConstraints< typename MatrixType::value_type > &constraints=AffineConstraints< typename MatrixType::value_type >())
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)
std::vector< unsigned int > serial(const std::vector< unsigned int > &targets, const std::function< RequestType(const unsigned int)> &create_request, const std::function< AnswerType(const unsigned int, const RequestType &)> &answer_request, const std::function< void(const unsigned int, const AnswerType &)> &process_answer, const MPI_Comm comm)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:470
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
*** braid_TestAll(app, comm_x, stdout, 0.0,(tstop-tstart)/ntime, *2 *(tstop-tstart)/ntime, my_Init, my_Free, my_Clone, *my_Sum, my_SpatialNorm, my_BufSize, my_BufPack, *my_BufUnpack, my_Coarsen, my_Interp, my_Residual, my_Step)
****code *  *  MPI_Finalize()
*braid_SplitCommworld & comm
void advance(std::tuple< I1, I2 > &t, const unsigned int n)