DiFfRG
|
In DiFfRG the computational setup is fully described by a numerical model. A model essentially is a description of a large set of coupled differential equations and some additional methods to handle data output.
In general, we have three components to any flow, which are FE functions \( u_i(x),\,i\in\{0,\dots,N_f-1\} \), variables \( v_a,\,a\in\{0,\dots,N_v-1\} \) and extractors \( e_b,\,b\in\{0,\dots,N_e-1\} \).
The latter two are just independent variables, whereas the FE functions depend additionaly on a field variable, \( u_i(x),\,x\in\mathbb{R}^d \). In other words, the FE functions explicitly live on a spatial discretization of the field space.
Any model you define should at least be inherited from the abstract model class DiFfRG::def::AbstractModel to ensure that all necessary methods are at least defined to do nothing, i.e.
Inside the class we can now overwrite all methods from DiFfRG::def::AbstractModel in order to implement the right system of flow equations.
The FE functions usually correspond to expansion coefficients in a derivative expansion. As an example consider a bosonic theory as in this paper: Treating a purely bosonic theory in first-order derivative expansion, the effective action is given by
\[\large \Gamma_k[\phi] = \int_x \bigg(\frac{1}{2}Z(\rho)(\partial_\mu\phi)^2 + V(\rho) \bigg)\,, \]
where \( \rho = \phi^2 / 2 \). A flowing reparametrization of the field \( \phi \) is being performed and is given by
\[\large \dot\phi(x) = \frac{1}{2} \eta(\rho) \phi\,. \]
where we introduced the anomalous dimension \( \eta = \frac{\partial_{t_+} Z}{Z} \).
The flow is then fully parametrized in terms of FE functions
\begin{align}\large u_1(x) &= m^2(\rho) = \partial_\rho V(\rho)\,, \\\large u_2(x) &= \eta(\rho)\,, \end{align}
where we also chose the field \( x = \rho \). We see here that the FE functions live on a spatial discretization of the d-dimensional field space \( \mathbb{R}^d \).
With the above ansatz one can quickly compute flow equations from the Wetterich equation,
\[\large k\partial_k \Gamma_k[\Phi] = \frac{1}{2}\text{Tr}\, G_{\alpha\beta}\,k\partial_k R^{\alpha\beta}\,. \]
We remark that the time
\[\large t = t_+ = \ln\left(\frac{\Lambda}{k}\right)\,,\]
as used in DiFfRG, is opposite in sign to the RG-time as defined in most literature, \(t_- = \ln\left(\frac{k}{\Lambda}\right)\). This is simply due to many time solvers not accepting negative time arguments.
In order to discretize the flow equations on a finite element space, the flow equations are expressed in the standard differential-algebraic form
\[\large m_i(\partial_t u_j, u_j, x) + \partial_x F_i(u_j, \partial_x u_j, \partial_x^2 u_j, e_b, v_a, x) + s_i(u_j, \partial_x u_j, \partial_x^2 u_j, e_b, v_a, x) = 0\,, \]
where \( m_i \) are called the mass functions, \( F_i \) the fluxes and \( s_i \) the sources. The latter two are functions of the FE functions, their derivatives, the field variable, the variables and the extractors.
In principle, the above system of equations can contain both equations containing the time derivatives, i.e. differential components, and equations without time derivatives, i.e. algebraic components. In order to solve the resulting DAEs one is currently restricted to the SUNDIALS IDA solver, which is however highly efficient and actually recommended for most cases.
Alternatively, the restricted formulation, allowing only for differential components,
\[\large m_{ij}(x) \partial_t u_j + \partial_x F_i(u_j, \partial_x u_j, \partial_x^2 u_j, e_b, v_a, x) + s_i(u_j, \partial_x u_j, \partial_x^2 u_j, e_b, v_a, x) = 0\,, \]
is used for all other provided ODE solvers, i.e. Runge-Kutta methods.
Note, that in the above definitions a change from \( t = t_+ \to t_- \) simply moves all terms onto the other side, i.e. when calculating the flow equations in the standard \(t_-\), one can still copy and paste everything without changing signs if the mass functions are simply \( m_i = \partial_{t_+} u_i = - \partial_{t_-} u_i \) (as is default).
The above components of standard form have direct analogues in the abstract numerical model which must be reimplemented for any system of flow equations. For actual implementation examples, especially regarding the template structure, please take a look at the models contained in the DiFfRG/models/
folder.
The relevant methods are also documented in DiFfRG::def::AbstractModel and read as follows:
m_i
argument is the resulting mass function \(m_i\), with \(N_f\) components. This method should fill the m_i
argument with the desired structure of the flow equation. x
is a d-dimensional array of field coordinates, and both u
(~ \(u_i(x)\)) and u_dot
(~ \(\partial_t u_i(x)\)) have \(N_f\) components. m_ij
argument is the resulting mass matrix \(m_{ij}\), with \(N_f\) components in each dimension. This method should fill the mass
argument with the desired structure of the flow equation. x
is a d-dimensional array of field coordinates. The standard implementation of this method simply sets \(m_{ij} = \delta_{ij}\).NumberType
and Solution
templates, whereas the rest can be dropped. F_i
has \(N_f\) components, x
gives the coordinate in field space and sol
contains all other arguments of the flux function. In practice, sol
is a std::tuple<...>
which containsNumberType
and Solution
templates, whereas the rest can be dropped. s_i
has \(N_f\) components, x
gives the coordinate in field space and sol
contains all other arguments of the flux function, with the layout as explained above in the flux case. The standard implementation of this method simply sets \(s_i = 0\).Picking up the example from above, we can now sketch the implementation of the numerical model as follows:
The actual numerical calculation of the flow equations (rather, their weak form) is done by the so-called assemblers. These are responsible for the actual discretization of the flow equations on the finite element space. In DiFfRG, we provide a set of assemblers for different discretizations, which are all derived from the abstract assembler class DiFfRG::AbstractAssembler.
Although all require at least the above interface methods, certain additional methods are required for certain discretizations. For example, the discontinuous Galerkin (DG) assemblers require the implementation of the numerical fluxes, and both discontinuous and continuous Galerkin (CG, also called simply FEM here) assemblers require the implementation of the boundary condition fluxes.
To understand the underlying numerics, see e.g. this review and also the excellent deal.ii tutorials.
For further reference, please refer to the documentation of the respective assemblers.
Underlying the assemblers are the actual discretizations, which are implemented in the DiFfRG::discretization namespace. These are responsible for the actual discretization of the field space, i.e. the construction of the finite element space. In DiFfRG, we provide a set of discretizations for different finite element spaces.
Putting everything together, we can write a straightforward main function to run the flow equations:
Here we have chosen a dDG
discretization, which is a discontinuous Galerkin discretization with a direct discontinuous Galerkin assembler. The timestepper is the SUNDIALS IDA solver, which is the recommended solver for most cases. For the dDG
discretization it is also necessary to supply a numerical flux, which can be done by modifying the numerical model as follows:
Here, the local Lax-Friedrichs flux has been used for the numerical fluxes and the boundaries have been defined to be inflow-/outflow. We have also chosen to use the autodiff functionality for the calculation of the jacobians.