318 * #include
"BraidFuncs.hh"
322 * This advances the solution forward by one time step.
323 * First some data is collected from the status
struct,
324 * namely the start and stop time and the current timestep
325 * number. The timestep size @f$\Delta t@f$ is calculated,
326 * and the step function from the HeatEquation is used to
330 *
int my_Step(braid_App app,
331 * braid_Vector ustop,
332 * braid_Vector fstop,
334 * braid_StepStatus status)
344 * braid_StepStatusGetLevel(status, &
level);
345 * braid_StepStatusGetTstartTstop(status, &tstart, &tstop);
346 * braid_StepStatusGetTIndex(status, &index);
348 * deltaT = tstop - tstart;
350 * ::Vector<double>& solution = u->data;
352 * HeatEquation<2>& heateq = app->eq;
354 * heateq.step(solution, deltaT, tstart, index);
362 * In
this function we initialize a vector at an arbitrary time.
363 * At
this point we don
't know anything about what the solution
364 * looks like, and we can really initialize to anything, so in
365 * this case use reinit to initialize the memory and set the
370 * my_Init(braid_App app,
372 * braid_Vector *u_ptr)
374 * my_Vector *u = new(my_Vector);
375 * int size = app->eq.size();
376 * u->data.reinit(size);
378 * app->eq.initialize(t, u->data);
387 * Here we need to copy the vector u into the vector v. We do this
388 * by allocating a new vector, then reinitializing the deal.ii
389 * vector to the correct size. The deal.ii reinitialization sets
390 * every value to zero, so next we need to iterate over the vector
391 * u and copy the values to the new vector v.
395 * my_Clone(braid_App app,
397 * braid_Vector *v_ptr)
400 * my_Vector *v = new(my_Vector);
401 * int size = u->data.size();
402 * v->data.reinit(size);
403 * for(size_t i=0, end=v->data.size(); i != end; ++i)
405 * v->data[i] = u->data[i];
414 * Here we need to free the memory used by vector u. This is
415 * pretty simple since the deal.ii vector is stored inside the
416 * XBraid vector, so we just delete the XBraid vector u and it
417 * puts the deal.ii vector out of scope and releases its memory.
421 * my_Free(braid_App app,
432 * This is to perform an axpy type operation. That is to say we
433 * do @f$y = \alpha x + \beta y@f$. Fortunately deal.ii already has
434 * this operation built in to its vector class, so we get the
435 * reference to the vector y and call the sadd method.
438 * int my_Sum(braid_App app,
445 * Vector<double>& vec = y->data;
446 * vec.sadd(beta, alpha, x->data);
453 * This calculates the spatial norm using the l2 norm. According
454 * to XBraid, this could be just about any spatial norm but we'll
455 * keep it simple and used deal.ii vector
's built in l2_norm method.
459 * my_SpatialNorm(braid_App app,
465 * dot = u->data.l2_norm();
473 * This function is called at various points depending on the access
474 * level specified when configuring the XBraid struct. This function
475 * is used to print out data during the run time, such as plots of the
476 * data. The status struct contains a ton of information about the
477 * simulation run. Here we get the current time and timestep number.
478 * The output_results function is called to plot the solution data.
479 * If the method of manufactured solutions is being used, then the
480 * error of this time step is computed and processed.
484 * my_Access(braid_App app,
486 * braid_AccessStatus astatus)
491 * braid_AccessStatusGetT(astatus, &t);
492 * braid_AccessStatusGetTIndex(astatus, &index);
494 * app->eq.output_results(index, t, u->data);
497 * if(index == app->final_step)
499 * app->eq.process_solution(t, index, u->data);
508 * This calculates the size of buffer needed to pack the solution
509 * data into a linear buffer for transfer to another processor via
510 * MPI. We query the size of the data from the HeatEquation class
511 * and return the buffer size.
515 * my_BufSize(braid_App app,
517 * braid_BufferStatus bstatus)
520 * int size = app->eq.size();
521 * *size_ptr = (size+1)*sizeof(double);
528 * This function packs a linear buffer with data so that the buffer
529 * may be sent to another processor via MPI. The buffer is cast to
530 * a type we can work with. The first element of the buffer is the
531 * size of the buffer. Then we iterate over solution vector u and
532 * fill the buffer with our solution data. Finally we tell XBraid
533 * how much data we wrote.
537 * my_BufPack(braid_App app,
540 * braid_BufferStatus bstatus)
544 * double *dbuffer = (double*)buffer;
545 * int size = u->data.size();
547 * for(int i=0; i != size; ++i)
549 * dbuffer[i+1] = (u->data)[i];
551 * braid_BufferStatusSetSize(bstatus, (size+1)*sizeof(double));
558 * This function unpacks a buffer that was received from a different
559 * processor via MPI. The size of the buffer is read from the first
560 * element, then we iterate over the size of the buffer and fill
561 * the values of solution vector u with the data in the buffer.
565 * my_BufUnpack(braid_App app,
567 * braid_Vector *u_ptr,
568 * braid_BufferStatus bstatus)
573 * my_Vector *u = NULL;
574 * double *dbuffer = (double*)buffer;
575 * int size = static_cast<int>(dbuffer[0]);
576 * u = new(my_Vector);
577 * u->data.reinit(size);
579 * for(int i = 0; i != size; ++i)
581 * (u->data)[i] = dbuffer[i+1];
590<a name="ann-src/BraidFuncs.hh"></a>
591<h1>Annotated version of src/BraidFuncs.hh</h1>
597 * /* -----------------------------------------------------------------------------
599 * * SPDX-License-Identifier: LGPL-2.1-or-later
600 * * Copyright (C) 2018 by Joshua Christopher
602 * * This file is part of the deal.II code gallery.
604 * * -----------------------------------------------------------------------------
607 * #ifndef _BRAIDFUNCS_H_
608 * #define _BRAIDFUNCS_H_
611 * * \file BraidFuncs.cc
612 * * \brief Contains the implementation of the mandatory X-Braid functions
614 * * X-Braid mandates several functions in order to drive the solution.
615 * * This file contains the implementation of said mandatory functions.
616 * * See the X-Braid documentation for more information.
617 * * There are several functions that are optional in X-Braid that may
618 * * or may not be implemented in here.
623 * /*-------- Third Party --------*/
624 * #include <deal.II/numerics/vector_tools.h>
627 * #include <braid_test.h>
629 * /*-------- Project --------*/
630 * #include "HeatEquation.hh"
634 * This struct contains all data that changes with time. For now
635 * this is just the solution data. When doing AMR this should
636 * probably include the triangulization, the sparsity pattern,
641 * * \brief Struct that contains the deal.ii vector.
643 * typedef struct _braid_Vector_struct
645 * ::Vector<double> data;
650 * This struct contains all the data that is unchanging with time.
654 * * \brief Struct that contains the HeatEquation and final
655 * * time step number.
657 * typedef struct _braid_App_struct
659 * HeatEquation<2> eq;
665 * * @brief my_Step - Takes a step in time, advancing the u vector
667 * * @param app - The braid app struct
668 * * @param ustop - The solution data at the end of this time step
669 * * @param fstop - RHS data (such as forcing function?)
670 * * @param u - The solution data at the beginning of this time step
671 * * @param status - Status structure that contains various info of this time
673 * * @return Success (0) or failure (1)
675 * int my_Step(braid_App app,
676 * braid_Vector ustop,
677 * braid_Vector fstop,
679 * braid_StepStatus status);
683 * * @brief my_Init - Initializes a solution data at the given time
684 * * For now, initializes the solution to zero no matter what time we are at
686 * * @param app - The braid app struct containing user data
687 * * @param t - Time at which the solution is initialized
688 * * @param u_ptr - The solution data that needs to be filled
690 * * @return Success (0) or failure (1)
693 * my_Init(braid_App app,
695 * braid_Vector *u_ptr);
699 * * @brief my_Clone - Clones a vector into a new vector
701 * * @param app - The braid app struct containing user data
702 * * @param u - The existing vector containing data
703 * * @param v_ptr - The empty vector that needs to be filled
705 * * @return Success (0) or failure (1)
708 * my_Clone(braid_App app,
710 * braid_Vector *v_ptr);
714 * * @brief my_Free - Deletes a vector
716 * * @param app - The braid app struct containing user data
717 * * @param u - The vector that needs to be deleted
719 * * @return Success (0) or failure (1)
722 * my_Free(braid_App app,
727 * * @brief my_Sum - Sums two vectors in an AXPY operation
728 * * The operation is y = alpha*x + beta*y
730 * * @param app - The braid app struct containing user data
731 * * @param alpha - The coefficient in front of x
732 * * @param x - A vector that is multiplied by alpha then added to y
733 * * @param beta - The coefficient of y
734 * * @param y - A vector that is multiplied by beta then summed with x
736 * * @return Success (0) or failure (1)
739 * my_Sum(braid_App app,
746 * * \brief Returns the spatial norm of the provided vector
748 * * Calculates and returns the spatial norm of the provided vector.
749 * * Interestingly enough, X-Braid does not specify a particular norm.
750 * * to keep things simple, we implement the Euclidean norm.
752 * * \param app - The braid app struct containing user data
753 * * \param u - The vector we need to take the norm of
754 * * \param norm_ptr - Pointer to the norm that was calculated, need to modify this
755 * * \return Success (0) or failure (1)
758 * my_SpatialNorm(braid_App app,
763 * * \brief Allows the user to output details
765 * * The Access function is called at various points to allow the user to output
766 * * information to the screen or to files.
767 * * The astatus parameter provides various information about the simulation,
768 * * see the XBraid documentation for details on what information you can get.
769 * * Example information is what the current timestep number and current time is.
770 * * If the access level (in parallel_in_time.cc) is set to 0, this function is
772 * * If the access level is set to 1, the function is called after the last
774 * * If the access level is set to 2, it is called every XBraid cycle.
776 * * \param app - The braid app struct containing user data
777 * * \param u - The vector containing the data at the status provided
778 * * \param astatus - The Braid status structure
779 * * \return Success (0) or failure (1)
782 * my_Access(braid_App app,
784 * braid_AccessStatus astatus);
787 * * \brief Calculates the size of a buffer for MPI data transfer
789 * * Calculates the size of the buffer that is needed to transfer
790 * * a solution vector to another processor.
791 * * The bstatus parameter provides various information on the
792 * * simulation, see the XBraid documentation for all possible
795 * * \param app - The braid app struct containing user data
796 * * \param size_ptr A pointer to the calculated size
797 * * \param bstatus The XBraid status structure
798 * * \return Success (0) or failure (1)
801 * my_BufSize(braid_App app,
803 * braid_BufferStatus bstatus);
806 * * \brief Linearizes a vector to be sent to another processor
808 * * Linearizes (packs) a data buffer with the contents of
809 * * some solution state u.
811 * * \param app - The braid app struct containing user data
812 * * \param u The vector that must be packed into buffer
813 * * \param buffer The buffer that must be filled with u
814 * * \param bstatus The XBraid status structure
815 * * \return Success (0) or failure (1)
818 * my_BufPack(braid_App app,
821 * braid_BufferStatus bstatus);
824 * * \brief Unpacks a vector that was sent from another processor
826 * * Unpacks a linear data buffer into the vector pointed to by
829 * * \param app - The braid app struct containing user data
830 * * \param buffer The buffer that must be unpacked
831 * * \param u_ptr The pointer to the vector that is filled
832 * * \param bstatus The XBraid status structure
833 * * \return Success (0) or failure (1)
836 * my_BufUnpack(braid_App app,
838 * braid_Vector *u_ptr,
839 * braid_BufferStatus bstatus);
841 * #endif // _BRAIDFUNCS_H_
845<a name="ann-src/HeatEquation.hh"></a>
846<h1>Annotated version of src/HeatEquation.hh</h1>
852 * /* -----------------------------------------------------------------------------
854 * * SPDX-License-Identifier: LGPL-2.1-or-later
855 * * Copyright (C) 2018 by Joshua Christopher
857 * * This file is part of the deal.II code gallery.
859 * * -----------------------------------------------------------------------------
862 * #ifndef _HEATEQUATION_H_
863 * #define _HEATEQUATION_H_
865 * #include <deal.II/base/utilities.h>
866 * #include <deal.II/base/quadrature_lib.h>
867 * #include <deal.II/base/function.h>
868 * #include <deal.II/base/logstream.h>
869 * #include <deal.II/lac/vector.h>
870 * #include <deal.II/lac/full_matrix.h>
871 * #include <deal.II/lac/dynamic_sparsity_pattern.h>
872 * #include <deal.II/lac/sparse_matrix.h>
873 * #include <deal.II/lac/solver_cg.h>
874 * #include <deal.II/lac/precondition.h>
875 * #include <deal.II/lac/affine_constraints.h>
876 * #include <deal.II/grid/tria.h>
877 * #include <deal.II/grid/grid_generator.h>
878 * #include <deal.II/grid/grid_refinement.h>
879 * #include <deal.II/grid/grid_out.h>
880 * #include <deal.II/grid/tria_accessor.h>
881 * #include <deal.II/grid/tria_iterator.h>
882 * #include <deal.II/dofs/dof_handler.h>
883 * #include <deal.II/dofs/dof_accessor.h>
884 * #include <deal.II/dofs/dof_tools.h>
885 * #include <deal.II/fe/fe_q.h>
886 * #include <deal.II/fe/fe_values.h>
887 * #include <deal.II/numerics/data_out.h>
888 * #include <deal.II/numerics/vector_tools.h>
889 * #include <deal.II/numerics/error_estimator.h>
890 * #include <deal.II/numerics/solution_transfer.h>
891 * #include <deal.II/numerics/matrix_tools.h>
892 * #include <deal.II/base/convergence_table.h>
896 * using namespace dealii;
900 * The HeatEquation class is describes the finite element
901 * solver for the heat equation. It contains all the functions
902 * needed to define the problem domain and advance the solution
912 * void step(Vector<double>& braid_data,
917 * int size() const; /// Returns the size of the solution vector
919 * void output_results(int a_time_idx,
921 * Vector<double>& a_solution) const;
923 * void initialize(double a_time,
924 * Vector<double>& a_vector) const;
926 * void process_solution(double a_time,
928 * const Vector<double>& a_vector);
931 * void setup_system();
932 * void solve_time_step(Vector<double>& a_solution);
934 * Triangulation<dim> triangulation;
936 * DoFHandler<dim> dof_handler;
938 * AffineConstraints<double> constraints;
940 * SparsityPattern sparsity_pattern;
941 * SparseMatrix<double> mass_matrix;
942 * SparseMatrix<double> laplace_matrix;
943 * SparseMatrix<double> system_matrix;
945 * Vector<double> system_rhs;
947 * std::ofstream myfile;
949 * const double theta;
953 * These were originally in the run() function but because
954 * I am splitting the run() function up into define and step
955 * they need to become member data
958 * Vector<double> tmp;
959 * Vector<double> forcing_terms;
961 * ConvergenceTable convergence_table;
966 * The RightHandSide class describes the RHS of the governing
967 * equations. In this case, it is the forcing function.
971 * class RightHandSide : public Function<dim>
980 * virtual double value (const Point<dim> &p,
981 * const unsigned int component = 0) const override;
984 * const double period;
989 * The BoundaryValues class describes the boundary conditions
990 * of the governing equations.
994 * class BoundaryValues : public Function<dim>
997 * virtual double value (const Point<dim> &p,
998 * const unsigned int component = 0) const override;
1003 * The RightHandSideMFG class describes the right hand side
1004 * function when doing the method of manufactured solutions.
1007 * template <int dim>
1008 * class RightHandSideMFG : public Function<dim>
1011 * virtual double value (const Point<dim> &p,
1012 * const unsigned int component = 0) const override;
1017 * The InitialValuesMFG class describes the initial values
1018 * when doing the method of manufactured solutions.
1021 * template <int dim>
1022 * class InitialValuesMFG : public Function<dim>
1025 * virtual double value (const Point<dim> &p,
1026 * const unsigned int component = 0) const override;
1031 * Provides the exact value for the manufactured solution. This
1032 * is used for the boundary conditions as well.
1035 * template <int dim>
1036 * class ExactValuesMFG : public Function<dim>
1040 * * \brief Computes the value at the given point and member data time
1042 * * Computes the exact value of the manufactured solution at point p and
1043 * * the member data time. See the class documentation and the design doc
1044 * * for details on what the exact solution is.
1046 * * \param p The point that the exact solution is computed at
1047 * * \param component The component of the exact solution (always 0 for now)
1048 * * \return double The exact value that was computed
1050 * virtual double value (const Point<dim> &p,
1051 * const unsigned int component = 0) const override;
1054 * * \brief Computes the gradient of the exact solution at the given point
1056 * * Computes the gradient of the exact/manufactured solution value at
1057 * * point p and member data time. See the design doc for details on
1058 * * what the gradient of the exact solution is
1060 * * \param p The point that the gradient is calculated at
1061 * * \param component The component of the system of equations this gradient is for
1062 * * \return Tensor<1,dim> A rank 1 tensor that contains the gradient
1063 * * in each spatial dimension
1065 * virtual Tensor<1,dim> gradient (const Point<dim> &p,
1066 * const unsigned int component = 0) const override;
1070 * #include "HeatEquationImplem.hh"
1072 * #endif // _HEATEQUATION_H_
1076<a name="ann-src/HeatEquationImplem.hh"></a>
1077<h1>Annotated version of src/HeatEquationImplem.hh</h1>
1083 * /* -----------------------------------------------------------------------------
1085 * * SPDX-License-Identifier: LGPL-2.1-or-later
1086 * * Copyright (C) 2018 by Joshua Christopher
1088 * * This file is part of the deal.II code gallery.
1090 * * -----------------------------------------------------------------------------
1093 * #include "Utilities.hh"
1095 * #include <iomanip>
1100 * Calculates the forcing function for the RightHandSide. See the
1101 * documentation for the math.
1104 * template <int dim>
1105 * double RightHandSide<dim>::value (const Point<dim> &p,
1106 * const unsigned int component) const
1109 * Assert (component == 0, ExcIndexRange(component, 0, 1));
1110 * Assert (dim == 2, ExcNotImplemented());
1112 * double time = this->get_time();
1114 * if ((p[0] > 0.5) && (p[1] > -0.5))
1116 * return std::exp(-0.5*(time-0.125)*(time-0.125)/(0.005));
1118 * else if ((p[0] > -0.5) && (p[1] > 0.5))
1120 * return std::exp(-0.5*(time-0.375)*(time-0.375)/(0.005));
1127 * return 0; // No forcing function
1132 * Calculates the forcing function for the method of manufactured
1133 * solutions. See the documentation for the math.
1136 * template <int dim>
1137 * double RightHandSideMFG<dim>::value (const Point<dim> &p,
1138 * const unsigned int component) const
1141 * Assert (component == 0, ExcIndexRange(component, 0, 1));
1142 * Assert (dim == 2, ExcNotImplemented());
1144 * double time = this->get_time();
1146 * double pi = numbers::PI;
1147 * return 4*pi*pi*std::exp(-4*pi*pi*time)*std::cos(2*pi*p[0])*std::cos(2*pi*p[1]);
1152 * Calculates the boundary conditions, essentially zero everywhere.
1155 * template <int dim>
1156 * double BoundaryValues<dim>::value (const Point<dim> &p,
1157 * const unsigned int component) const
1161 * Assert (component == 0, ExcIndexRange(component, 0, 1));
1167 * Calculates the exact solution (and thus also boundary conditions)
1168 * for the method of manufactured solutions.
1171 * template <int dim>
1172 * double ExactValuesMFG<dim>::value (const Point<dim> &p,
1173 * const unsigned int component) const
1176 * Assert (component == 0, ExcIndexRange(component, 0, 1));
1178 * double time = this->get_time();
1179 * const double pi = numbers::PI;
1181 * return std::exp(-4*pi*pi*time)*std::cos(2*pi*p[0])*std::cos(2*pi*p[1]);
1186 * Calculates the gradient of the exact solution for the method of manufactured
1187 * solutions. See the documentation for the math.
1190 * template <int dim>
1191 * Tensor<1,dim> ExactValuesMFG<dim>::gradient (const Point<dim> &p,
1192 * const unsigned int) const
1194 * Assert (dim == 2, ExcNotImplemented());
1196 * Tensor<1,dim> return_value;
1197 * const double pi = numbers::PI;
1198 * double time = this->get_time();
1199 * return_value[0] = -2*pi*std::exp(-4*pi*pi*time)*std::cos(2*pi*p[1])*std::sin(2*pi*p[0]);
1200 * return_value[1] = -2*pi*std::exp(-4*pi*pi*time)*std::cos(2*pi*p[0])*std::sin(2*pi*p[1]);
1201 * return return_value;
1206 * Calculates the initial values for the method of manufactured solutions.
1207 * See the documentation for the math.
1210 * template <int dim>
1211 * double InitialValuesMFG<dim>::value (const Point<dim> &p,
1212 * const unsigned int component) const
1215 * Assert (component == 0, ExcIndexRange(component, 0, 1));
1216 * const double pi = numbers::PI;
1218 * return std::cos(2*pi*p[0])*std::cos(2*pi*p[1]);
1221 * template <int dim>
1222 * HeatEquation<dim>::HeatEquation ()
1225 * dof_handler(triangulation),
1230 * template <int dim>
1231 * void HeatEquation<dim>::initialize(double a_time,
1232 * Vector<double>& a_vector) const
1237 * We only initialize values in the manufactured solution case
1240 * InitialValuesMFG<dim> iv_function;
1241 * iv_function.set_time(a_time);
1242 * VectorTools::project (dof_handler, constraints,
1243 * QGauss<dim>(fe.degree+1), iv_function,
1251 * If not the MFG solution case, a_vector is already zero'd so
do nothing
1256 *
template <
int dim>
1257 *
void HeatEquation<dim>::setup_system()
1259 * dof_handler.distribute_dofs(fe);
1261 * constraints.clear ();
1264 * constraints.close();
1271 * sparsity_pattern.copy_from(dsp);
1274 * laplace_matrix.reinit(sparsity_pattern);
1275 * system_matrix.reinit(sparsity_pattern);
1284 * system_rhs.reinit(dof_handler.n_dofs());
1288 *
template <
int dim>
1289 *
void HeatEquation<dim>::solve_time_step(
Vector<double>& a_solution)
1291 *
SolverControl solver_control(1000, 1e-8 * system_rhs.l2_norm());
1295 * preconditioner.
initialize(system_matrix, 1.0);
1297 * cg.solve(system_matrix, a_solution, system_rhs,
1300 * constraints.distribute(a_solution);
1305 *
template <
int dim>
1306 *
void HeatEquation<dim>::output_results(
int a_time_idx,
1312 * vtk_flags.
time = a_time;
1313 * vtk_flags.cycle = a_time_idx;
1318 * data_out.attach_dof_handler(dof_handler);
1319 * data_out.add_data_vector(a_solution,
"U");
1321 * data_out.build_patches();
1323 *
const std::string filename =
"solution-"
1326 * std::ofstream output(filename.c_str());
1327 * data_out.write_vtk(output);
1332 * We define the geometry here,
this is called on each processor
1333 * and doesn
't change in time. Once doing AMR, this won't need
1337 *
template <
int dim>
1338 *
void HeatEquation<dim>::define()
1340 *
const unsigned int initial_global_refinement = 6;
1347 * tmp.reinit (dof_handler.n_dofs());
1348 * forcing_terms.reinit (dof_handler.n_dofs());
1353 * Here we
advance the solution forward in time. This is done
1354 * the same way as in the
loop in @ref step_26
"step-26"'s run function.
1358 * void HeatEquation<dim>::step(Vector<double>& braid_data,
1366 * mass_matrix.vmult(system_rhs, braid_data);
1368 * laplace_matrix.vmult(tmp, braid_data);
1370 * system_rhs.add(-(1 - theta) * deltaT, tmp);
1373 * RightHandSideMFG<dim> rhs_function;
1375 * RightHandSide<dim> rhs_function;
1377 * rhs_function.set_time(a_time);
1378 * VectorTools::create_right_hand_side(dof_handler,
1379 * QGauss<dim>(fe.degree+1),
1383 * forcing_terms = tmp;
1384 * forcing_terms *= deltaT * theta;
1386 * rhs_function.set_time(a_time - deltaT);
1387 * VectorTools::create_right_hand_side(dof_handler,
1388 * QGauss<dim>(fe.degree+1),
1392 * forcing_terms.add(deltaT * (1 - theta), tmp);
1393 * system_rhs += forcing_terms;
1395 * system_matrix.copy_from(mass_matrix);
1396 * system_matrix.add(theta * deltaT, laplace_matrix);
1398 * constraints.condense (system_matrix, system_rhs);
1404 * If we are doing the method of manufactured solutions
1405 * then we set the boundary conditions to the exact solution.
1406 * Otherwise the boundary conditions are zero.
1409 * ExactValuesMFG<dim> boundary_values_function;
1411 * BoundaryValues<dim> boundary_values_function;
1413 * boundary_values_function.set_time(a_time);
1415 * std::map<types::global_dof_index, double> boundary_values;
1416 * VectorTools::interpolate_boundary_values(dof_handler,
1418 * boundary_values_function,
1421 * MatrixTools::apply_boundary_values(boundary_values,
1427 * solve_time_step(braid_data);
1431 * int HeatEquation<dim>::size() const
1433 * return dof_handler.n_dofs();
1438 * This function computes the error for the time step when doing
1439 * the method of manufactured solutions. First the exact values
1440 * is calculated, then the difference per cell is computed for
1441 * the various norms, and the error is computed. This is written
1442 * out to a pretty table.
1445 * template<int dim> void
1446 * HeatEquation<dim>::process_solution(double a_time,
1448 * const Vector<double>& a_vector)
1452 * Compute the exact value for the manufactured solution case
1455 * ExactValuesMFG<dim> exact_function;
1456 * exact_function.set_time(a_time);
1458 * Vector<double> difference_per_cell (triangulation.n_active_cells());
1459 * VectorTools::integrate_difference(dof_handler,
1462 * difference_per_cell,
1463 * QGauss<dim>(fe.degree+1),
1464 * VectorTools::L2_norm);
1466 * const double L2_error = VectorTools::compute_global_error(triangulation,
1467 * difference_per_cell,
1468 * VectorTools::L2_norm);
1470 * VectorTools::integrate_difference(dof_handler,
1473 * difference_per_cell,
1474 * QGauss<dim>(fe.degree+1),
1475 * VectorTools::H1_seminorm);
1477 * const double H1_error = VectorTools::compute_global_error(triangulation,
1478 * difference_per_cell,
1479 * VectorTools::H1_seminorm);
1481 * const QTrapezoid<1> q_trapez;
1482 * const QIterated<dim> q_iterated (q_trapez, 5);
1483 * VectorTools::integrate_difference (dof_handler,
1486 * difference_per_cell,
1488 * VectorTools::Linfty_norm);
1489 * const double Linfty_error = VectorTools::compute_global_error(triangulation,
1490 * difference_per_cell,
1491 * VectorTools::Linfty_norm);
1493 * const unsigned int n_active_cells = triangulation.n_active_cells();
1494 * const unsigned int n_dofs = dof_handler.n_dofs();
1496 * pout() << "Cycle " << a_index << ':
'
1498 * << " Number of active cells: "
1501 * << " Number of degrees of freedom: "
1505 * convergence_table.add_value("cycle", a_index);
1506 * convergence_table.add_value("cells", n_active_cells);
1507 * convergence_table.add_value("dofs", n_dofs);
1508 * convergence_table.add_value("L2", L2_error);
1509 * convergence_table.add_value("H1", H1_error);
1510 * convergence_table.add_value("Linfty", Linfty_error);
1512 * convergence_table.set_precision("L2", 3);
1513 * convergence_table.set_precision("H1", 3);
1514 * convergence_table.set_precision("Linfty", 3);
1516 * convergence_table.set_scientific("L2", true);
1517 * convergence_table.set_scientific("H1", true);
1518 * convergence_table.set_scientific("Linfty", true);
1520 * convergence_table.set_tex_caption("cells", "\\# cells");
1521 * convergence_table.set_tex_caption("dofs", "\\# dofs");
1522 * convergence_table.set_tex_caption("L2", "@fL^2@f-error");
1523 * convergence_table.set_tex_caption("H1", "@fH^1@f-error");
1524 * convergence_table.set_tex_caption("Linfty", "@fL^\\infty@f-error");
1526 * convergence_table.set_tex_format("cells", "r");
1527 * convergence_table.set_tex_format("dofs", "r");
1529 * std::cout << std::endl;
1530 * convergence_table.write_text(std::cout);
1532 * std::ofstream error_table_file("tex-conv-table.tex");
1533 * convergence_table.write_tex(error_table_file);
1538<a name="ann-src/Utilities.cc"></a>
1539<h1>Annotated version of src/Utilities.cc</h1>
1545 * /* -----------------------------------------------------------------------------
1547 * * SPDX-License-Identifier: LGPL-2.1-or-later
1548 * * Copyright (C) 2018 by Joshua Christopher
1550 * * This file is part of the deal.II code gallery.
1552 * * -----------------------------------------------------------------------------
1555 * #include "Utilities.hh"
1558 * #include <fstream>
1567 * The shared variables
1573 * static std::string s_pout_filename ;
1574 * static std::string s_pout_basename ;
1575 * static std::ofstream s_pout ;
1577 * static bool s_pout_init = false ;
1578 * static bool s_pout_open = false ;
1582 * in parallel, compute the filename give the basename
1583 * [NOTE: dont call this before MPI is initialized.]
1586 * static void setFileName()
1588 * static const size_t ProcnumSize = 1 + 10 + 1 ; //'.
' + 10digits + '\0
'
1589 * char procnum[ProcnumSize] ;
1590 * snprintf( procnum ,ProcnumSize ,".%d" ,procID);
1591 * s_pout_filename = s_pout_basename + procnum ;
1596 * in parallel, close the file if nec., open it and check for success
1599 * static void openFile()
1601 * if ( s_pout_open )
1605 * s_pout.open( s_pout_filename.c_str() );
1608 * if open() fails, we have problems, but it's better
1609 * to
try again later than to make believe it succeeded
1612 * s_pout_open = (
bool)s_pout ;
1616 * std::ostream& pout()
1621 * the common
case is _open ==
true, which just returns s_pout
1624 *
if ( ! s_pout_open )
1628 * the uncommon cae: the file isn
't opened, MPI may not be
1629 * initialized, and the basename may not have been set
1632 * int flag_i, flag_f;
1633 * MPI_Initialized(&flag_i);
1634 * MPI_Finalized(&flag_f);
1637 * app hasn't
set a basename yet, so
set the
default
1640 *
if ( ! s_pout_init )
1642 * s_pout_basename =
"pout" ;
1643 * s_pout_init = true ;
1647 *
if MPI not initialized, we can
't open the file so return cout
1650 * if ( ! flag_i || flag_f)
1652 * return std::cout; // MPI hasn't been started yet, or has ended....
1656 *
MPI is initialized, so file must not be, so open it
1663 *
finally, in
case the open failed,
return cout
1666 *
if ( ! s_pout_open )
1668 *
return std::cout ;
1679<a name=
"ann-src/Utilities.hh"></a>
1680<h1>Annotated version of src/
Utilities.hh</h1>
1696 * #ifndef _UTILITIES_H_
1697 * #define _UTILITIES_H_
1699 * #include <iostream>
1703 * This preprocessor macro is used on function arguments
1704 * that are not used in the function. It is used to
1705 * suppress compiler warnings.
1708 * #define UNUSED(x) (void)(x)
1712 * Contains the current
MPI processor ID.
1715 *
extern int procID;
1719 *
Function to
return the ostream to write out to. In
MPI
1720 * mode it returns a stream to a file named pout.<#> where
1721 * <#> is the procID. This allows the user to write output
1722 * from each processor to a separate file. In
serial mode
1723 * (no
MPI), it returns the standard output.
1726 * std::ostream& pout();
1731<a name=
"ann-src/parallel_in_time.cc"></a>
1732<h1>Annotated version of src/parallel_in_time.cc</h1>
1750 * #include
"BraidFuncs.hh"
1751 * #include
"HeatEquation.hh"
1752 * #include
"Utilities.hh"
1754 * #include <fstream>
1755 * #include <iostream>
1757 *
int main(
int argc,
char *argv[])
1761 *
using namespace dealii;
1766 * MPI_Init(&argc, &argv);
1767 *
comm = MPI_COMM_WORLD;
1768 * MPI_Comm_rank(
comm, &rank);
1778 *
double tstart = 0.0;
1779 *
double tstop = 0.002;
1781 * my_App *app =
new(my_App);
1783 * braid_Init(MPI_COMM_WORLD,
comm, tstart, tstop, ntime, app,
1784 * my_Step, my_Init, my_Clone, my_Free, my_Sum, my_SpatialNorm,
1785 * my_Access, my_BufSize, my_BufPack, my_BufUnpack, &core);
1789 *
int max_levels = 3;
1796 *
double tol = 1.e-7;
1805 *
int min_coarse = 10;
1809 *
int wrapper_tests = 0;
1812 *
int print_level = 1;
1813 *
int access_level = 1;
1814 *
int use_sequential= 0;
1816 * braid_SetPrintLevel( core, print_level);
1817 * braid_SetAccessLevel( core, access_level);
1818 * braid_SetMaxLevels(core, max_levels);
1821 * braid_SetMinCoarse( core, min_coarse );
1822 * braid_SetSkip(core, skip);
1823 * braid_SetNRelax(core, -1, nrelax);
1826 * braid_SetAbsTol(core, tol);
1829 * braid_SetCFactor(core, -1, cfactor);
1832 * braid_SetMaxIter(core, max_iter);
1833 * braid_SetSeqSoln(core, use_sequential);
1836 * app->final_step = ntime;
1838 * braid_Drive(core);
1842 * Free the memory now that we are done
1845 * braid_Destroy(core);
1852 * MPI_Comm_free(&
comm);
1857 *
catch (std::exception &exc)
1859 * std::cerr << std::endl << std::endl
1860 * <<
"----------------------------------------------------"
1862 * std::cerr <<
"Exception on processing: " << std::endl << exc.what()
1863 * << std::endl <<
"Aborting!" << std::endl
1864 * <<
"----------------------------------------------------"
1871 * std::cerr << std::endl << std::endl
1872 * <<
"----------------------------------------------------"
1874 * std::cerr <<
"Unknown exception!" << std::endl <<
"Aborting!"
1876 * <<
"----------------------------------------------------"
1887<a name=
"ann-test/test_braid.cc"></a>
1888<h1>Annotated version of test/test_braid.cc</h1>
1904 * #include
"BraidFuncs.hh"
1906 * #include <braid.h>
1907 * #include <braid_test.h>
1909 * #include <iostream>
1911 *
int main(
int argc,
char** argv)
1915 * MPI_Init(&argc, &argv);
1916 *
comm = MPI_COMM_WORLD;
1917 * MPI_Comm_rank(
comm, &rank);
1919 * my_App *app =
new(my_App);
1922 *
double time = 0.2;
1924 * braid_Int init_access_result = braid_TestInitAccess(app,
1931 * (
void)init_access_result;
1933 * braid_Int clone_result = braid_TestClone(app,
1941 * (
void)clone_result;
1943 * braid_Int sum_result = braid_TestSum(app,
1954 * braid_Int norm_result = braid_TestSpatialNorm(app,
1963 * (
void)norm_result;
1965 * braid_Int buf_result = braid_TestBuf(app,
1981 * braid_SplitCommworld(&
comm, 1, &comm_x, &comm_t);
1985 * 2*(tstop-tstart)/ntime, my_Init, my_Free, my_Clone,
1986 * my_Sum, my_SpatialNorm, my_BufSize, my_BufPack,
1987 * my_BufUnpack, my_Coarsen, my_Interp, my_Residual, my_Step);
void set_flags(const FlagType &flags)
void initialize(const MatrixType &A, const AdditionalData ¶meters=AdditionalData())
void refine_global(const unsigned int times=1)
__global__ void set(Number *val, const Number s, const size_type N)
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())
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)
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.)
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)
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)
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)