From 15b07a6c9f291dffe71b40e9cf56a214869e1308 Mon Sep 17 00:00:00 2001 From: Joshua Christopher Date: Fri, 18 May 2018 22:03:41 -0600 Subject: [PATCH] Added a ton more documentation, proofread the Readme, removed some extra cmake lines that weren't needed --- parallel_in_time/CMakeLists.txt | 23 -- parallel_in_time/Readme.md | 236 ++++++++++++--------- parallel_in_time/src/BraidFuncs.cc | 59 +++++- parallel_in_time/src/BraidFuncs.hh | 13 +- parallel_in_time/src/HeatEquation.hh | 24 ++- parallel_in_time/src/HeatEquationImplem.hh | 74 ++++--- parallel_in_time/src/Utilities.cc | 2 +- parallel_in_time/src/Utilities.hh | 9 + 8 files changed, 252 insertions(+), 188 deletions(-) diff --git a/parallel_in_time/CMakeLists.txt b/parallel_in_time/CMakeLists.txt index 5474895..db5df0f 100644 --- a/parallel_in_time/CMakeLists.txt +++ b/parallel_in_time/CMakeLists.txt @@ -60,29 +60,6 @@ endif() # Include the braid paths and libraries include_directories(${BRAID_DIR}) - -############################## -# Finally start building stuff - -# First, build the "library" that consists of all the source files -# set(CMAKE_LIBRARY_OUTPUT_DIRECTORY ${LIB_PATH}) -# set(CMAKE_ARCHIVE_OUTPUT_DIRECTORY ${LIB_PATH}) -# add_library(${LIB_NAME} ${SRC}) - - -# Next build the main function -# set(CMAKE_RUNTIME_OUTPUT_DIRECTORY ${BIN_PATH}) -# add_executable(${BIN_NAME} ${MAIN_SRC}) -# target_link_libraries(${MAIN_NAME} ${LIB_NAME}) # link my library -# target_link_libraries(${MAIN_NAME} ${BRAID_DIR}libbraid.a) # link braid -# DEAL_II_SETUP_TARGET(${MAIN_NAME}) # Link dealii -# if(USE_MPI) -# set_target_properties(${MAIN_NAME} PROPERTIES -# LINK_FLAGS "${MPI_LINK_FLAGS}") # Link MPI -# set_target_properties(${MAIN_NAME} PROPERTIES -# COMPILE_FLAGS "${MPI_COMPILE_FLAGS}") # Use MPI compile flags -# endif(USE_MPI) - DEAL_II_INITIALIZE_CACHED_VARIABLES() PROJECT(${TARGET}) DEAL_II_INVOKE_AUTOPILOT() diff --git a/parallel_in_time/Readme.md b/parallel_in_time/Readme.md index bdc5f7d..db5af9b 100644 --- a/parallel_in_time/Readme.md +++ b/parallel_in_time/Readme.md @@ -1,6 +1,5 @@ -# Parallel in Time deal.ii - -## Overview +Overview +======== Over the last few years, the clock speed per processor core has stagnated. This stagnation has lead to the design of larger core counts in high performance computing machines. @@ -9,65 +8,69 @@ Perhaps the largest bottleneck in concurrency for time-dependent simulations is Traditional time integration methods solve for the time domain sequentially, and as the spatial grid is refined a proportionally larger number of time steps must be taken to maintain accuracy and stability constraints. While solving the time domain sequentially with a traditional time integration method is an optimal algorithm of order $\mathcal{O}(n)$, the $n$ time steps are not solved concurrently. -The goal of this project is to make use of the XBraid library from Lawrence Livermore National Laboratory to solve the time domain in parallel using multigrid-reduction-in-time techniques. +The goal of this project is to make use of the XBraid library from Lawrence Livermore National Laboratory to solve the time domain in parallel using multigrid reduction in time techniques. The XBraid library is implemented in C and aims to be a non-intrusive method to implement parallel time marching methods into existing codes. -## Implementation +Implementation +============== -### XBraid introduction +XBraid introduction +------------------- In order to use the XBraid library, several data structures and functions must be implemented and provided to the XBraid solver struct. -The two required data structures are the app and vector structures. -In general, the app struct contains the time independent data and the vector struct contains the time dependent data. +The two required data structures are the app and vector structures. +In general, the app struct contains the time independent data and the vector struct contains the time dependent data. For this initial example, the time independent data includes the mesh which is fixed for all time steps, and the time dependent data is the solution state vector. -The functions tell XBraid how to perform operations on the data type used by your solver, in this case deal.ii uses the Vector data type. +The functions tell XBraid how to perform operations on the data type used by your solver, in this case deal.ii uses the Vector data type. These operations include how to initialize the data at a given time, how to sum the data, and how to pack and unpack linear buffers for transmission to other processors via MPI. The XBraid documentation should be read for a full list of functions that must be implemented and the details of what the function should do. -The typical format is the function is called with arguments of the app struct, one or more vector structs, and a status struct that contains information on the current status of the XBraid simulation (the current multigrid iteration, the level the function is being called from, the time and timestep number, etc.). -Perhaps the most important function is the step function. -This function tells XBraid how to advance the solution forward in time from the initial to the final times given in the status struct. +The typical format is the function is called with arguments of the app struct, one or more vector structs, and a status struct that contains information on the current status of the XBraid simulation (the current multigrid iteration, the level the function is being called from, the time and timestep number, etc.). + +Perhaps the most important function is the step function. +This function tells XBraid how to advance the solution forward in time from the initial to the final times given in the status struct. This method uses a traditional time integration method such as the fourth order explicit Runge Kutta method. -### deal.ii details +deal.ii details +--------------- -The solver used in this example is based off the heat equation solver from the [Step-26 Tutorial](https://dealii.org/developer/doxygen/deal.II/step_26.html . -The HeatEquation class becomes member data to XBraid's app struct, and XBraid's vector struct becomes a wrapper for deal.ii's Vector data type. -The HeatEquation class cannot simply be used as is though as it contains both time dependent and time independent member data. +The solver used in this example is based off the heat equation solver from the step-26 tutorial of deal.ii. +The HeatEquation class becomes member data to XBraid’s app struct, and XBraid’s vector struct becomes a wrapper for deal.ii’s Vector data type. +The HeatEquation class cannot simply be used as is though as it contains both time dependent and time independent member data. In order to simplify the problem the adaptive mesh refinement is removed. Theoretically XBraid is capable of working with adaptive mesh refinement and in fact contains support for time refinement (which is also not used for simplicity). -All adaptive mesh refinement functionality is removed from the sovler. -The time-dependent solution state vectors are also removed from the HeatEquation member data. -That time-dependent data will be provided at each timestep by XBraid via the vector struct. +All adaptive mesh refinement functionality is removed from the solver. +The time-dependent solution state vectors are also removed from the HeatEquation member data. +That time-dependent data will be provided at each timestep by XBraid via the vector struct. -### Math details - -In the default mode, this code solves the heat equation with, +Governing Equations {#math-details} +------------------- +In the default mode, this code solves the heat equation, @f{align} -\frac{\partial u}{\partial t} - \Delta u = f(\boldsymbol{x},t), \qquad \forall\boldsymbol{x}\in\Omega,t\in\left( 0,T \right), + \frac{\partial u}{\partial t} - \Delta u = f(\boldsymbol{x},t), \qquad \forall\boldsymbol{x}\in\Omega,t\in\left( + 0,T \right), @f} with initial conditions, @f{align} -u(\boldsymbol{x},0) = u_0(\boldsymbol{x}) = 0, \qquad \forall \boldsymbol{x}\in\Omega, + u(\boldsymbol{x},0) = u_0(\boldsymbol{x}) = 0, \qquad \forall \boldsymbol{x}\in\Omega, @f} -and dirichlet boundary conditions, +and Dirichlet boundary conditions, @f{align} -u(\boldsymbol{x},t) = g(\boldsymbol{x},t) = 0, \qquad \forall \boldsymbol{x}\in\partial\Omega,t\in\left( 0,T \right) + u(\boldsymbol{x},t) = g(\boldsymbol{x},t) = 0, \qquad \forall \boldsymbol{x}\in\partial\Omega,t\in\left( 0,T \right) @f} and forcing function, @f{align} - f(\mathbf x, t) = - \left\{ - \begin{array}{ll} - \chi_1(t) & \text{if \(x>0.5\) and \(y>-0.5\)} \\ - \chi_2(t) & \text{if \(x>-0.5\) and \(y>0.5\)} - \end{array} - \right. + f(\mathbf x, t) = \left\{ + \begin{array}{ll} + \chi_1(t) & \text{if \(x>0.5\) and \(y>-0.5\)} \\ + \chi_2(t) & \text{if \(x>-0.5\) and \(y>0.5\)} + \end{array} + \right. @f} with, @f{align} - \chi_1(t) = \exp\left(-0.5\frac{(t-0.125)^2}{0.005}\right) + \chi_1(t) = \exp\left(-0.5\frac{(t-0.125)^2}{0.005}\right) @f} and, @f{align} @@ -75,69 +78,82 @@ and, @f} for some time $t$. The forcing function is a Gaussian pulse in time that is centered around 0.125 time units for $\chi_1$ and 0.375 time units for $\chi_2$. -A gaussian function was chosen because it is a continuous function in time and so the solution state can be compared bitwise with the serial-in-time solution. +A Gaussian function was chosen because it is a continuous function in time and so the solution state can be compared bit-wise with the serial-in-time solution. + +Method of Manufactured Solutions +-------------------------------- The method of manufactured solutions is used to test the correctness of the implementation. In the method of manufactured solutions, we create a solution $u_h$ to the heat equation, then compute the boundary conditions, initial conditions, and forcing functions required to generate that solution. +This method is explained further in the step-7 tutorial of deal.ii. The created solution used is, @f{align} -u_h = \exp\left(-4\pi^2t\right)\cos(2 \pi x)\cos(2 \pi y), \qquad \forall \boldsymbol{x} \in \Omega \cup \partial\Omega + u_h = \exp\left( -4\pi^2t \right) \cos(2 \pi x) \cos(2 \pi y), \qquad \forall \boldsymbol{x} \in \Omega \cup \partial\Omega @f} with derivatives, @f{align} -\frac{\partial u}{\partial t} &= -4 \pi^2 e^{-4\pi^2t}\cos(2 \pi x)\cos(2 \pi y), \\ --\Delta u &= 8 \pi^2 \exp\left(-4\pi^2t\right)\cos(2 \pi x)\cos(2 \pi y) \\ -\frac{\partial u}{\partial x} &= -2 \pi \exp\left(-4\pi^2t\right)\sin(2\pi x)\cos(2\pi y) \\ -\frac{\partial u}{\partial x} &= -2 \pi \exp\left(-4\pi^2t\right)\cos(2\pi x)\sin(2\pi y) + \frac{\partial u}{\partial t} &= -4 \pi^2 \exp{-4\pi 2t} \cos(2 \pi x)\cos(2 \pi y), \\ + -\Delta u &= 8 \pi^2 \exp\left(-4\pi^2t\right)\cos(2 \pi x)\cos(2 \pi y) \\ + \frac{\partial u}{\partial x} &= -2\pi \exp\left(-4\pi^2t\right)\sin(2\pi x)\cos(2\pi y) \\ + \frac{\partial u}{\partial x} &= -2 \pi \exp\left(-4\pi^2t\right)\cos(2\pi x)\sin(2\pi y) @f} and therefore we specify the forcing term, initial conditions, and boundary conditions of the governing equations as, @f{align} -f(\boldsymbol{x},t) &= 4 \pi^2 \exp\left(-4\pi^2t\right)\cos(2 \pi x)\cos(2 \pi y), &\forall\boldsymbol{x}\in\Omega,t\in\left( 0,T \right), \\ -u_0(\boldsymbol{x}) &= \cos(2 \pi x)\cos(2 \pi y), &\forall \boldsymbol{x}\in\Omega, \\ -g(\boldsymbol{x},t) &= \exp\left(-4\pi^2t\right)\cos(2 \pi x)\cos(2 \pi y), &\forall \boldsymbol{x} \in \partial\Omega. + f(\boldsymbol{x},t) &= 4 \pi^2 \exp\left(-4\pi^2t\right)\cos(2 \pi x)\cos(2 \pi y), &&\forall\boldsymbol{x}\in\Omega,t\in\left( 0,T \right), \\ + u_0(\boldsymbol{x}) &= \cos(2 \pi x)\cos(2\pi y), &&\forall \boldsymbol{x}\in\Omega, \\ + g(\boldsymbol{x},t) &= \exp\left(-4\pi^2t\right)\cos(2 \pi x)\cos(2\pi y), &&\forall \boldsymbol{x} \in \partial\Omega. @f} The manufactured solution is run on progressively more refined grids and the solution generated by the finite element method is compared to the exact solution $u_h$. The convergence rate of the error is calculated with @f{align} - \Delta \epsilon_n = \frac{\ln{\epsilon_{n-1}/\epsilon_{n}}}{\ln{r_n}} + \Delta \epsilon_n = \frac{\ln{\epsilon_{n-1}/\epsilon_{n}}}{\ln{\frac\Delta x_{n-1}/\Delta x_{n}}} @f} -where $\Delta \epsilon_n$ is the convergence rate of error $\epsilon$ between a mesh $n$ and -coarser mesh $n-1$ that have a refinement ratio of $r_n$. -Shown in Table 1 is the convergence rate as the mesh is refined. +where $\Delta \epsilon_n$ is the convergence rate of error $\epsilon$ between a mesh $n$ and coarser mesh $n-1$ that have a refinement ratio of $r_n$. +Shown in Table 1 are the errors of the $\textrm{L}^2$, $\textrm{H}^1$, and $\textrm{L}^\infty$ norms as the mesh is refined, and shown in Table 2 is the convergence rate. The $\Delta t$ is reduced by a factor of 2 for every global refinement of the mesh. @f{align} -\begin{tabular}{|c|c|c|c|c|c|} \hline - cycle & \# cells & \# dofs & L^2-error & H^1-error & L^\infty-error \\ \hline - 125 & 48 & 65 & 6.036e-03 & 6.970e-02 & 7.557e-03\\ \hline - 250 & 192 & 225 & 1.735e-03 & 3.414e-02 & 2.721e-03 \\ \hline - 500 & 768 & 833 & 4.513e-04 & 1.690e-02 & 7.410e-04 \\ \hline - 1000 & 3072 & 3201 & 1.140e-04 & 8.426e-03 & 1.877e-04 \\ \hline - 2000 & 12288 & 12545 & 2.859e-05 & 4.209e-03 & 4.715e-05 \\ \hline -\end{tabular} + \begin{tabular}{|c|c|c|c|c|c|} \hline + cycles & \# cells & \# dofs & $\textrm{L}^2$-error & $\textrm{H}^1$-error & $\textrm{L}^\infty$-error \\ \hline + 125 & 48 & 65 & 6.036e-03 & 6.970e-02 & 7.557e-03\\ \hline + 250 & 192 & 225 & 1.735e-03 & 3.414e-02 & 2.721e-03 \\ \hline + 500 & 768 & 833 & 4.513e-04 & 1.690e-02 & 7.410e-04 \\ \hline + 1000 & 3072 & 3201 & 1.140e-04 & 8.426e-03 & 1.877e-04 \\ \hline + 2000 & 12288 & 12545 & 2.859e-05 & 4.209e-03 & 4.715e-05 \\ \hline + \end{tabular} @f} @f{align} -\begin{tabular}{|c|c|c|c|c|c|} \hline - cycle & \# cells & \# dofs & Slope L^2 & Slope H^1 & Slope \textrm{L}^\infty \\ \hline - 125 & 48 & 65 & --- & --- & --- \\ \hline - 250 & 192 & 225 & 1.798 & 1.030 & 1.474 \\ \hline - 500 & 768 & 833 & 1.943 & 1.014 & 1.877 \\ \hline - 1000 & 3072 & 3201 & 1.985 & 1.004 & 1.981 \\ \hline - 2000 & 12288 & 12545 & 1.995 & 1.001 & 1.993 \\ \hline -\end{tabular} + \begin{tabular}{|c|c|c|c|c|c|} \hline + cycles & \# cells & \# dofs & Slope $\textrm{L}^2$ & Slope $\textrm{H}^1$ & Slope $\textrm{L}^\infty$ \\ \hline + 125 & 48 & 65 & --- & --- & --- \\ \hline + 250 & 192 & 225 & 1.798 & 1.030 & 1.474 \\ \hline + 500 & 768 & 833 & 1.943 & 1.014 & 1.877 \\ \hline + 1000 & 3072 & 3201 & 1.985 & 1.004 & 1.981 \\ \hline + 2000 & 12288 & 12545 & 1.995 & 1.001 & 1.993 \\ \hline + \end{tabular} @f} + + +The convergence rate is plotted as a function of the grid spacing in the figure below. As can be seen, the slope converges at a second order rate for the $\textrm{L}^2$ and $\textrm{L}^\infty$ norms and at a first order rate for the $\textrm{H}^1$ norms. This is expected behavior as second order finite elements are used. -### Code Organization +![Grid convergence results](./doc/results/gridConvergence.png) + +Code Organization +----------------- #### The src directory -The entry point of the code is in parallel_in_time.cc and sets up XBraid for a simulation. The XBraid setup involves initializing the app struct and configuring XBraid for the desired number of timesteps, number of iterations, etc. -The functions implemented for XBraid's use are declared in BraidFuncs.hh and defined in BraidFuncs.cc. The HeatEquation class and all deal.ii functionality is declared in HeatEquation.hh and defiend in HeatEquationImplem.hh. Since HeatEquation is a class template, its definition file HeatEquationImplem.hh is included at the bottom of HeatEquation.hh. Lastly various helper functions and variables such as the current processor id and the output stream are declared in Utilities.hh and defined in Utilities.cc. +The entry point of the code is in parallel\_in\_time.cc and sets up XBraid for a simulation. +The XBraid setup involves initializing the app struct and configuring XBraid for the desired number of timesteps, number of iterations, and so forth. +The functions implemented for XBraid’s use are declared in `BraidFuncs.hh` and defined in `BraidFuncs.cc`. +The HeatEquation class and all deal.ii functionality is declared in `HeatEquation.hh` and defined in `HeatEquationImplem.hh`. +Since HeatEquation is a class template, its definition file `HeatEquationImplem.hh` is included at the bottom of `HeatEquation.hh`. +Lastly various helper functions and variables such as the current processor id and the output stream are declared in `Utilities.hh` and defined in `Utilities.cc`. #### The test directory @@ -147,77 +163,89 @@ These tests verify the correct implementation of the various functions. #### The doc directory This directory is for storing further documentation of the code. -Not much is in this directory right now as most of the documentation is in this file (Readme.md) or in comments in the source code files. +Not much is in this directory right now as most of the documentation is in the Readme or in comments of the source code files. +Documentation is generated from the Readme and code comments by deal.ii’s documentation process. -## Compiling +Compiling +========= -To compile, you need deal.ii and XBraid to be installed with development headers somwehere on your system. +To compile, you need deal.ii and XBraid to be installed with development headers somewhere on your system. Some implementation of MPI such as OpenMPI with development headers must also be installed. -The source code for deal.ii is available at [https://dealii.org/](https://dealii.org/) and the source code for XBraid is available at [https://computation.llnl.gov/projects/parallel-time-integration-multigrid](https://computation.llnl.gov/projects/parallel-time-integration-multigrid). +The source code for deal.ii is available at [deal.ii’s website](https://dealii.org/) and the source code for XBraid is available at [LLNL’s website](https://computation.llnl.gov/projects/parallel-time-integration-multigrid). See the documentation of each package for compilation and installation instructions. -Depending on where they are installed, parallel_in_time may need help finding these libraries. -To find deal.ii, parallel_in_time first looks in typical deal.ii install directories followed by one directory up(../), two directories up (../../), and lastly in the environment variable DEAL_II_DIR. -In contrast, XBraid currently does not have any default locations to look for and so the environment variable BRAID_DIR must be specified. -For MPI, parallel_in_time looks in standard installation folders only, for that reason I recommend you install MPI with your package manager. +Depending on where they are installed, `parallel_in_time` may need help finding these libraries. +To find deal.ii, `parallel_in_time` first looks in typical deal.ii install directories followed by one directory up (`../`), two directories up (`../../`), and lastly in the environment variable `DEAL_II_DIR`. +In contrast, XBraid currently does not have any default locations to look for and so the environment variable `BRAID_DIR` must be specified. +For MPI, `parallel_in_time` looks in standard installation folders only, for that reason I recommend you install MPI with your package manager. A compile process of the parallel in time code may look like, - mkdir build - cd build - BRAID_DIR=/path/to/braid/ cmake ../ - make +``` + mkdir build + cd build + BRAID_DIR=/path/to/braid/ cmake ../ + make +``` -There is currently no option to install parallel_in_time anywhere. -The binaries are generated in the bin folder, and tests are placed into the test folder. -Options that can be passed to CMake for parallel_in_time include: +There is currently no option to install parallel\_in\_time anywhere. +The binaries are generated in the bin folder, and tests are placed into the test folder. +Options that can be passed to CMake for parallel\_in\_time include: -* CMAKE_BUILD_TYPE=Debug/Release -* DO_MFG=ON/OFF -* USE_MPI=ON/OFF +``` + CMAKE_BUILD_TYPE=Debug/Release + DO_MFG=ON/OFF + USE_MPI=ON/OFF +``` -The build type specifies whether to compile with debugging symbols, assertions, and optimizations or not. +The build type specifies whether to compile with debugging symbols, +assertions, and optimizations or not. -The option for manufactured solutions (DO_MFG) switches from solving the "standard" heat equation to solving a known solution heat equation so that the correctness of the code can be tested. +The option for manufactured solutions (DO\_MFG) switches from solving the “standard” heat equation to solving a known solution heat equation so that the correctness of the code can be tested. -Lastly the MPI option is only used to specify where to write output information when using the pout() function from the Utilities.hh file. -If USE_MPI is set to ON, then every processor writes to its own file called pout.<#> where <#> is the processor number. -If USE_MPI is set to OFF, then every processor writes to stdout. +Lastly the MPI option is only used to specify where to write output information when using the pout() function from the `Utilities.hh` file. +If USE\_MPI is set to ON, then every processor writes to its own file called pout.<\#> where <\#> is the processor number. +If USE\_MPI is set to OFF, then every processor writes to stdout. -## Running +Running +======= -Once parallel_in_time has been compiled, the program can be run by calling the binary generated in ./build/bin/. -The test can be run by calling ctest from inside the build/ directory. -Unless the output path has been changed in the source code (currently hardcoded), then the output files will be placed into the folder the command was called from. +Once parallel\_in\_time has been compiled, the program can be run by calling the binary generated in `./build/bin/`. +The test can be run by calling ctest from inside the `build/` directory. +Unless the output path has been changed in the source code (currently hard coded), then the output files will be placed into the folder the command was called from. There are no argument parameters for any of the binaries or tests. -## Results +Results +======= -To test the performance, a strong scaling study is run. +To test the performance, a strong scaling study is run. The spatial grid is fixed at 3201 degrees of freedom, and the spatial grid consists of 25,000 uniform timesteps. No spatial parallelization is used and the grid is fixed for all timesteps. -The parallel in time solution is solved using XBraid's multigrid reduction in time algorithm on 1, 2, 4, 16, 32, and 64 processors. +The parallel in time solution is solved using XBraid’s multigrid reduction in time algorithm on 1, 2, 4, 16, 32, and 64 processors. The serial in time solution is run on a single processor using traditional sequential time stepping. -The results are shown in the figure below. -Running the multigrid algorithm on a single processor takes about an order of magnitude longer to run on a single processor than the serial algorihtm on a single processor. +The results are shown in Figure \[fig:strongscaling\]. +Running the multigrid algorithm on a single processor takes about an order of magnitude longer to run on a single processor than the serial algorithm on a single processor. At 16 processors the multigrid algorithm the wall clock time is approximately the same for the serial algorithm as for the multigrid algorithm, and for 32 and 64 processors in time the wall clock is faster by about a factor of 2 for 64 processors. ![Strong scaling results](./doc/results/strongscaling.png) -## Conclusions +Conclusions +=========== -Using 64 times as many proceossors results in a speedup factor of approximately 2. +Using 64 times as many processors results in a speedup factor of approximately 2. This is a fairly heavy computational cost for a slight speedup. -For comparison, in the reference paper Falgout et. al. they achieved a speedup of approximately 10 times when using 256 times as many processors as their serial case when solving the heat equation in two dimensions. +For comparison, in the reference paper they achieved a speedup of approximately 10 times when using 256 times as many processors as their serial case when solving the heat equation in two dimensions. A similar increase in processors may be needed to increase the speedup factor of this code from 2 to 10. -The choice of whether to use serial in time or parallel in time largely rests in the size of problem being solved and the amount of computing resources available. Increasing the required number of timesteps will benefit the parallel in time algorithm provided enough extra resources are available. +The choice of whether to use serial in time or parallel in time largely rests in the size of problem being solved, the amount of computing resources available, and the urgency of the problem being solved. +Increasing the required number of timesteps will benefit the parallel in time algorithm provided enough extra resources are available. -## Future Work +Future Work +=========== There are many routes available for future work. First there are many optimizations that could be made to speed up the existing code base such as exploring the different multigrid cycles and finding the optimal load balancing. -Ultimately those optimizations will probably only result in marginal gains per the Falgout paper. +Ultimately those optimizations will probably only result in marginal gains per the reference paper. Allowing XBraid to prolong and restrict the spatial grid may be one of the more promising avenues of improvement. Future work that is of interest to the authors of XBraid is the development of adaptive mesh refinement (AMR) in a parallel in time algorithm. diff --git a/parallel_in_time/src/BraidFuncs.cc b/parallel_in_time/src/BraidFuncs.cc index 2160fbc..92cd5e2 100644 --- a/parallel_in_time/src/BraidFuncs.cc +++ b/parallel_in_time/src/BraidFuncs.cc @@ -1,6 +1,11 @@ -/*-------- Project --------*/ #include "BraidFuncs.hh" +// This advances the solution forward by one time step. +// First some data is collected from the status struct, +// namely the start and stop time and the current timestep +// number. The timestep size $\Delta t$ is calculated, +// and the step function from the HeatEquation is used to +// advance the solution. int my_Step(braid_App app, braid_Vector ustop, braid_Vector fstop, @@ -30,13 +35,12 @@ int my_Step(braid_App app, return 0; } -/** - * In this function we initialize a vector at an arbitrary time. - * At this point we don't know anything about what the solution - * looks like, and we can really initialize to anything, so in - * this case use reinit to initialize the memory and set the - * values to zero. - */ + +// In this function we initialize a vector at an arbitrary time. +// At this point we don't know anything about what the solution +// looks like, and we can really initialize to anything, so in +// this case use reinit to initialize the memory and set the +// values to zero. int my_Init(braid_App app, double t, @@ -53,6 +57,11 @@ my_Init(braid_App app, return 0; } +// Here we need to copy the vector u into the vector v. We do this +// by allocating a new vector, then reinitializing the deal.ii +// vector to the correct size. The deal.ii reinitialization sets +// every value to zero, so next we need to iterate over the vector +// u and copy the values to the new vector v. int my_Clone(braid_App app, braid_Vector u, @@ -71,6 +80,10 @@ my_Clone(braid_App app, return 0; } +// Here we need to free the memory used by vector u. This is +// pretty simple since the deal.ii vector is stored inside the +// XBraid vector, so we just delete the XBraid vector u and it +// puts the deal.ii vector out of scope and releases its memory. int my_Free(braid_App app, braid_Vector u) @@ -81,6 +94,10 @@ my_Free(braid_App app, return 0; } +// This is to perform an axpy type operation. That is to say we +// do $y = \alpha x + \beta y$. Fortunately deal.ii already has +// this operation built in to its vector class, so we get the +// reference to the vector y and call the sadd method. int my_Sum(braid_App app, double alpha, braid_Vector x, @@ -94,6 +111,9 @@ int my_Sum(braid_App app, return 0; } +// This calculates the spatial norm using the l2 norm. According +// to XBraid, this could be just about any spatial norm but we'll +// keep it simple and used deal.ii vector's built in l2_norm method. int my_SpatialNorm(braid_App app, braid_Vector u, @@ -107,6 +127,14 @@ my_SpatialNorm(braid_App app, return 0; } +// This function is called at various points depending on the access +// level specified when configuring the XBraid struct. This function +// is used to print out data during the run time, such as plots of the +// data. The status struct contains a ton of information about the +// simulation run. Here we get the current time and timestep number. +// The output_results function is called to plot the solution data. +// If the method of manufactured solutions is being used, then the +// error of this time step is computed and processed. int my_Access(braid_App app, braid_Vector u, @@ -123,7 +151,6 @@ my_Access(braid_App app, #if DO_MFG if(index == app->final_step) { - pout() << "Doing error calc of step: " << index << std::endl; app->eq.process_solution(t, index, u->data); } #endif @@ -131,6 +158,10 @@ my_Access(braid_App app, return 0; } +// This calculates the size of buffer needed to pack the solution +// data into a linear buffer for transfer to another processor via +// MPI. We query the size of the data from the HeatEquation class +// and return the buffer size. int my_BufSize(braid_App app, int *size_ptr, @@ -143,6 +174,12 @@ my_BufSize(braid_App app, return 0; } +// This function packs a linear buffer with data so that the buffer +// may be sent to another processor via MPI. The buffer is cast to +// a type we can work with. The first element of the buffer is the +// size of the buffer. Then we iterate over soltuion vector u and +// fill the buffer with our solution data. Finally we tell XBraid +// how much data we wrote. int my_BufPack(braid_App app, braid_Vector u, @@ -163,6 +200,10 @@ my_BufPack(braid_App app, return 0; } +// This function unpacks a buffer that was recieved from a different +// processor via MPI. The size of the buffer is read from the first +// element, then we iterate over the size of the buffer and fill +// the values of solution vector u with the data in the buffer. int my_BufUnpack(braid_App app, void *buffer, diff --git a/parallel_in_time/src/BraidFuncs.hh b/parallel_in_time/src/BraidFuncs.hh index 0032d89..9e84c21 100644 --- a/parallel_in_time/src/BraidFuncs.hh +++ b/parallel_in_time/src/BraidFuncs.hh @@ -13,7 +13,6 @@ * */ -/*-------- System --------*/ /*-------- Third Party --------*/ #include @@ -24,15 +23,23 @@ /*-------- Project --------*/ #include "HeatEquation.hh" +// This struct contains all data that changes with time. For now +// this is just the solution data. When doing AMR this should +// probably include the triangulization, the sparsity patter, +// constraints, etc. /** - * \brief Struct that contains the data + * \brief Struct that contains the deal.ii vector. */ typedef struct _braid_Vector_struct { dealii::Vector data; } my_Vector; -// Wrap the Heat Equation in a struct +// This struct contains all the data that is unchanging with time. +/** + * \brief Struct that contains the HeatEquation and final + * time step number. + */ typedef struct _braid_App_struct { HeatEquation<2> eq; diff --git a/parallel_in_time/src/HeatEquation.hh b/parallel_in_time/src/HeatEquation.hh index 07737e4..b743c3a 100644 --- a/parallel_in_time/src/HeatEquation.hh +++ b/parallel_in_time/src/HeatEquation.hh @@ -34,6 +34,10 @@ using namespace dealii; +// The HeatEquation class is describes the finite element +// solver for the heat equation. It contains all the functions +// needed to define the problem domain and advance the solution +// in time. template class HeatEquation { @@ -47,9 +51,6 @@ public: int size() const; /// Returns the size of the solution vector - void dump_vec(std::ofstream& file, - const Vector& vector) const; - void output_results(int a_time_idx, double a_time, Vector& a_solution) const; @@ -91,6 +92,8 @@ private: ConvergenceTable convergence_table; }; +// The RightHandSide class describes the RHS of the governing +// equations. In this case, it is the forcing function. template class RightHandSide : public Function { @@ -108,6 +111,8 @@ private: const double period; }; +// The BoundaryValues class describes the boundary conditions +// of the governing equations. template class BoundaryValues : public Function { @@ -116,6 +121,8 @@ public: const unsigned int component = 0) const; }; +// The RightHandSideMFG class describes the right hand side +// function when doing the method of manufactured solutions. template class RightHandSideMFG : public Function { @@ -124,7 +131,8 @@ public: const unsigned int component = 0) const; }; - +// The InitialValuesMFG class describes the initial values +// when doing the method of manufactured solutions. template class InitialValuesMFG : public Function { @@ -133,12 +141,8 @@ public: const unsigned int component = 0) const; }; -/** - * Provides the exact value for the manufactured solution. - * User must set the time using ExactValuesMFG::set_time(const Number new_time); - * Calculates \f$u_h=e^(-4\pi\pi t)*\cos(2\pi x)\cos(2\pi y)\f$ for all points given to it. - * In the context of manufactured solutions, this is used on the boundary. - */ +// Provides the exact value for the manufactured solution. This +// is used for the boundary conditions as well. template class ExactValuesMFG : public Function { diff --git a/parallel_in_time/src/HeatEquationImplem.hh b/parallel_in_time/src/HeatEquationImplem.hh index 050fec6..3b18f11 100644 --- a/parallel_in_time/src/HeatEquationImplem.hh +++ b/parallel_in_time/src/HeatEquationImplem.hh @@ -2,6 +2,8 @@ #include #include "Utilities.hh" +// Calculates the forcing function for the RightHandSide. See the +// documentation for the math. template double RightHandSide::value (const Point &p, const unsigned int component) const @@ -28,6 +30,8 @@ double RightHandSide::value (const Point &p, return 0; // No forcing function } +// Calculates the forcing function for the method of manufactured +// solutions. See the documentation for the math. template double RightHandSideMFG::value (const Point &p, const unsigned int component) const @@ -38,11 +42,11 @@ double RightHandSideMFG::value (const Point &p, double time = this->get_time(); - // return the manufactured solution of the right hand side double pi = numbers::PI; return 4*pi*pi*std::exp(-4*pi*pi*time)*std::cos(2*pi*p[0])*std::cos(2*pi*p[1]); } +// Calculates the boundary conditions, essentially zero everywhere. template double BoundaryValues::value (const Point &p, const unsigned int component) const @@ -53,6 +57,8 @@ double BoundaryValues::value (const Point &p, return 0; } +// Calculates the exact solution (and thus also boundary conditions) +// for the method of manufactured solutions. template double ExactValuesMFG::value (const Point &p, const unsigned int component) const @@ -62,16 +68,18 @@ double ExactValuesMFG::value (const Point &p, double time = this->get_time(); const double pi = numbers::PI; - // Return our manufactured solution boundary value + return std::exp(-4*pi*pi*time)*std::cos(2*pi*p[0])*std::cos(2*pi*p[1]); } - -// TODO: Find bug in here +// Calculates the gradient of the exact solution for the method of manufactured +// solutions. See the documentation for the math. template Tensor<1,dim> ExactValuesMFG::gradient (const Point &p, const unsigned int) const { + Assert (dim == 2, ExcNotImplemented()); + Tensor<1,dim> return_value; const double pi = numbers::PI; double time = this->get_time(); @@ -80,6 +88,8 @@ Tensor<1,dim> ExactValuesMFG::gradient (const Point &p, return return_value; } +// Calculates the initial values for the method of manufactured solutions. +// See the documentation for the math. template double InitialValuesMFG::value (const Point &p, const unsigned int component) const @@ -87,7 +97,7 @@ double InitialValuesMFG::value (const Point &p, (void) component; Assert (component == 0, ExcIndexRange(component, 0, 1)); const double pi = numbers::PI; - // Return our manufactured solution initial value + return std::cos(2*pi*p[0])*std::cos(2*pi*p[1]); } @@ -123,15 +133,6 @@ void HeatEquation::setup_system() { dof_handler.distribute_dofs(fe); - // pout() << std::endl - // << "===========================================" - // << std::endl - // << "Number of active cells: " << triangulation.n_active_cells() - // << std::endl - // << "Number of degrees of freedom: " << dof_handler.n_dofs() - // << std::endl - // << std::endl; - constraints.clear (); DoFTools::make_hanging_node_constraints (dof_handler, constraints); @@ -201,7 +202,9 @@ void HeatEquation::output_results(int a_time_idx, data_out.write_vtk(output); } -// This function won't make much sense in real parallel in time... +// We define the geometry here, this is called on each processor +// and doesn't change in time. Once doing AMR, this won't need +// to exist anymore. template void HeatEquation::define() { @@ -216,14 +219,14 @@ void HeatEquation::define() forcing_terms.reinit (dof_handler.n_dofs()); } +// Here we advance the solution forward in time. This is done +// the same way as in the loop in step-26's run function. template void HeatEquation::step(Vector& braid_data, double deltaT, double a_time, int a_time_idx) { - // Set old solution to the braid data - // old_solution = braid_data; a_time += deltaT; ++a_time_idx; @@ -263,7 +266,9 @@ void HeatEquation::step(Vector& braid_data, { #if DO_MFG - // Set boundary to exact value in MFG solution + // If we are doing the method of manufactured solutions + // then we set the boundary conditions to the exact solution. + // Otherwise the boundary conditions are zero. ExactValuesMFG boundary_values_function; #else BoundaryValues boundary_values_function; @@ -291,18 +296,11 @@ int HeatEquation::size() const return dof_handler.n_dofs(); } -template void -HeatEquation::dump_vec(std::ofstream& file, - const Vector& vector) const -{ - file << "Dumping vec:"; - for(int i=0; i != vector.size(); ++i) - { - file << "\n" << vector[i]; - } - file << std::endl; -} - +// This function computes the error for the time step when doing +// the method of manufactured solutions. First the exact values +// is calculated, then the difference per cell is computed for +// the various norms, and the error is computed. This is written +// out to a pretty table. template void HeatEquation::process_solution(double a_time, int a_index, @@ -350,14 +348,14 @@ HeatEquation::process_solution(double a_time, const unsigned int n_active_cells = triangulation.n_active_cells(); const unsigned int n_dofs = dof_handler.n_dofs(); - std::cout << "Cycle " << a_index << ':' - << std::endl - << " Number of active cells: " - << n_active_cells - << std::endl - << " Number of degrees of freedom: " - << n_dofs - << std::endl; + pout() << "Cycle " << a_index << ':' + << std::endl + << " Number of active cells: " + << n_active_cells + << std::endl + << " Number of degrees of freedom: " + << n_dofs + << std::endl; convergence_table.add_value("cycle", a_index); convergence_table.add_value("cells", n_active_cells); diff --git a/parallel_in_time/src/Utilities.cc b/parallel_in_time/src/Utilities.cc index e0019eb..dca1c3d 100644 --- a/parallel_in_time/src/Utilities.cc +++ b/parallel_in_time/src/Utilities.cc @@ -7,7 +7,7 @@ int procID = 0; -// the shared variables +// The shared variables static std::string s_pout_filename ; static std::string s_pout_basename ; diff --git a/parallel_in_time/src/Utilities.hh b/parallel_in_time/src/Utilities.hh index 6e41cd0..607fa47 100644 --- a/parallel_in_time/src/Utilities.hh +++ b/parallel_in_time/src/Utilities.hh @@ -2,9 +2,18 @@ #define _UTILITIES_H_ #include +// This preprocessor macro is used on function arguments +// that are not used in the function. It is used to +// suppress compiler warnings. #define UNUSED(x) (void)(x) +// Contains the current MPI processor ID. extern int procID; +// Function to return the ostream to write out to. In MPI +// mode it returns a stream to a file named pout.<#> where +// <#> is the procID. This allows the user to write output +// from each processor to a separate file. In serial mode +// (no MPI), it returns the standard output. std::ostream& pout(); #endif -- 2.39.5