Go back to the Contents page.


Press Show to reveal the code chunks.


# Create a clipboard button on the rendered HTML page
source(here::here("clipboard.R")); clipboard
# Set seed for reproducibility
set.seed(1982) 
# Set global options for all code chunks
knitr::opts_chunk$set(
  # Disable messages printed by R code chunks
  message = FALSE,    
  # Disable warnings printed by R code chunks
  warning = FALSE,    
  # Show R code within code chunks in output
  echo = TRUE,        
  # Include both R code and its results in output
  include = TRUE,     
  # Evaluate R code chunks
  eval = TRUE,       
  # Enable caching of R code chunks for faster rendering
  cache = FALSE,      
  # Align figures in the center of the output
  fig.align = "center",
  # Enable retina display for high-resolution figures
  retina = 2,
  # Show errors in the output instead of stopping rendering
  error = TRUE,
  # Do not collapse code and output into a single block
  collapse = FALSE
)
# Start the figure counter
fig_count <- 0
# Define the captioner function
captioner <- function(caption) {
  fig_count <<- fig_count + 1
  paste0("Figure ", fig_count, ": ", caption)
}
library(MetricGraph)
library(ggplot2)
library(reshape2)
library(plotly)
capture.output(
  knitr::purl(here::here("functionality.Rmd"), output = here::here("functionality.R")),
  file = here::here("purl_log.txt")
)
source(here::here("functionality.R"))

1 Optimal control of fractional diffusion equations on metric graphs

1.1 Problem statement

Let \(\Gamma = (\mathcal{V},\mathcal{E})\) be a metric graph. Let \(u_d: \Gamma \times(0, T) \rightarrow \mathbb{R}\) be the desired state and \(\mu>0\) a regularization parameter. We define the cost functional

\[\begin{equation} \label{eq:costfun} \tag{1} J(u, z)=\frac{1}{2} \int_0^T\left(\left\|u-u_d\right\|_{L_2(\Gamma)}^2+\mu\|z\|_{L_2(\Gamma)}^2\right) dt \end{equation}\]

Let \(f:\Gamma\times (0,T)\rightarrow\mathbb{R}\) and \(u_0: \Gamma \rightarrow \mathbb{R}\) be fixed functions. We will call them right-hand side and initial datum, respectively. Let \(\alpha\in(0,2]\) and \(z: \Gamma \times(0, T) \rightarrow \mathbb{R}\) denote the control variable. We shall be concerned with the following PDE-constrained optimization problem: Find

\[\begin{equation} \label{eq:min_pro} \tag{2} \min\; J(u, z) \end{equation}\] subject to the fractional diffusion equation \[\begin{equation} \label{eq:maineq} \tag{3} \left\{ \begin{aligned} \partial_t u(s,t) + (\kappa^2 - \Delta_\Gamma)^{\alpha/2} u(s,t) &= f(s,t)+z(s,t), && \quad (s,t) \in \Gamma \times (0, T), \\ u(s,0) &= u_0(s), && \quad s \in \Gamma, \end{aligned} \right. \end{equation}\] with \(u(\cdot,t)\) satisfying the Kirchhoff vertex conditions \[\begin{equation} \label{eq:Kcond} \tag{4} \mathcal{K} = \left\{\phi\in C(\Gamma)\;\middle|\; \forall v\in \mathcal{V}:\; \sum_{e\in\mathcal{E}_v}\partial_e \phi(v)=0 \right\} \end{equation}\] and the control constraints \[\begin{align} \label{control_constraints} \tag{5} a(s,t)\leq z(s,t)\leq b(s,t)\;\mathrm{a.e.}\;(s,t)\in\Gamma \times(0, T). \end{align}\]

1.2 Numerical scheme

The optimal variables \((\bar{u}, \bar{p}, \bar{z})\) satisfy

\[\begin{equation} \label{eq:maineqoptimal} \tag{6} \left\{ \begin{aligned} \partial_t \bar{u}(s,t) + (\kappa^2 - \Delta_\Gamma)^{\alpha/2} \bar{u}(s,t) &= f(s,t)+\bar{z}(s,t), && \quad (s,t) \in \Gamma \times (0, T), \\ \bar{u}(s,0) &= u_0(s), && \quad s \in \Gamma, \end{aligned} \right. \end{equation}\] and \[\begin{equation} \label{eq:adjointeq} \tag{7} \left\{ \begin{aligned} -\partial_t \bar{p}(s,t) + (\kappa^2 - \Delta_\Gamma)^{\alpha/2} \bar{p}(s,t) &= \bar{u}(s,t)-u_d(s,t), && \quad (s,t) \in \Gamma \times (0, T), \\ \bar{p}(s,T) &= 0, && \quad s \in \Gamma, \end{aligned} \right. \end{equation}\] with \[\begin{align} \label{zz} \tag{8} \bar{z}(s,t) = \max\left\{a(s,t),\min\left\{b(s,t),-\dfrac{1}{\mu}\bar{p}(s,t)\right\}\right\}. \end{align}\] By considering the change of variable \(\bar{q}(s,t) = \bar{p}(s,T-t)\), the fractional adjoint problem \(\eqref{eq:adjointeq}\) becomes a forward-in-time problem where the transformed adjoint state \(\bar{q}\) solves \[\begin{equation} \label{transformed_adjoint_state} \tag{9} \left\{ \begin{aligned} \partial_t \bar{q}(s,t) + (\kappa^2 - \Delta_\Gamma)^{\alpha/2} \bar{q}(s,t) &= \bar{v}(s,t)-v_d(s,t), && \quad (s,t) \in \Gamma \times (0, T), \\ \bar{q}(s,0) &= 0, && \quad s \in \Gamma, \end{aligned} \right. \end{equation}\] since \(\partial_t\bar{q}(s,t) = -\partial_t\bar{p}(s,T-t)\) and \(\bar{q}(s, 0) =\bar{p}(s,T-0)= \bar{p}(s,T)=0\). Here, \(\bar{v}(s,t) = \bar{u}(s,T-t)\) and \(v_d(s,t) = u_d(s,T-t)\).

Given \(\bar{u}\) and \(u_d\), we can time-reverse them to obtain \(\bar{v}\) and \(v_d\) and then use the same numerical scheme we use for the forward problem \(\eqref{eq:maineqoptimal}\) to solve the adjoint problem \(\eqref{transformed_adjoint_state}\). The control variable \(\bar{z}\) is then computed using \(\eqref{zz}\).

The numerical scheme for \(\eqref{eq:maineqoptimal}\) and \(\eqref{transformed_adjoint_state}\) are given by (see the Functionality page) \[\begin{align} \label{numericalscheme1} \tag{10} \bar{\boldsymbol{U}}^{k+1} = \boldsymbol{C}^{-1}\left(\sum_{k=1}^{m+1} a_k\left(\dfrac{\boldsymbol{L}\boldsymbol{C}^{-1}}{\kappa^2}-p_k\boldsymbol{I}\right)^{-1} + r\boldsymbol{I}\right) (\boldsymbol{C}\bar{\boldsymbol{U}}^k+\tau (\boldsymbol{F}^{k+1}+\bar{\boldsymbol{Z}}^{k+1})), \end{align}\] and \[\begin{align} \label{thenumericalscheme2} \tag{11} \bar{\boldsymbol{Q}}^{k+1} = \boldsymbol{C}^{-1}\left(\sum_{k=1}^{m+1} a_k\left(\dfrac{\boldsymbol{L}\boldsymbol{C}^{-1}}{\kappa^2}-p_k\boldsymbol{I}\right)^{-1} + r\boldsymbol{I}\right) (\boldsymbol{C}\bar{\boldsymbol{Q}}^k+\tau (\bar{\boldsymbol{V}}^{k+1}-\boldsymbol{V_d}^{k+1})), \end{align}\] with initial conditions \(\bar{\boldsymbol{U}}^0 = \boldsymbol{U}_0\) and \(\bar{\boldsymbol{Q}}^0 = \boldsymbol{0}\). Observe that \(\eqref{numericalscheme1}\)-\(\eqref{thenumericalscheme2}\) is a coupled problem.

For notation convenience, let

  • \(\bar{\boldsymbol{U}}\) be the matrix with columns \(\bar{\boldsymbol{U}}^k\) obtained from the numerical scheme \(\eqref{numericalscheme1}\)
  • \(\boldsymbol{F}\) be the matrix with columns \(\boldsymbol{F}^k\) and entries \(\boldsymbol{F}_{ik}=(\psi^i_h,f(\cdot,t_k))_{L_2(\Gamma)}\)
  • \(\bar{\boldsymbol{Z}}\) be the matrix with columns \(\bar{\boldsymbol{Z}}^k\) and entries \(\bar{\boldsymbol{Z}}_{ik}=(\psi^i_h,\bar{z}(\cdot,t_k))_{L_2(\Gamma)}\)
  • \(\bar{\boldsymbol{Q}}\) be the matrix with columns \(\bar{\boldsymbol{Q}}^k\) obtained from the numerical scheme \(\eqref{thenumericalscheme2}\)
  • \(\bar{\boldsymbol{V}}\) be the matrix with columns \(\bar{\boldsymbol{V}}^k\) and entries \(\bar{\boldsymbol{V}}^{ik}=(\psi^i_h,\bar{v}(\cdot,t_k))_{L_2(\Gamma)} =(\psi^i_h,\bar{u}(\cdot,T-t_k))_{L_2(\Gamma)}\)
  • \(\boldsymbol{V_d}\) be the matrix with columns \(\boldsymbol{V_d}^k\) and entries \(\boldsymbol{V_d}^{ik}=(\psi^i_h,v_d(\cdot,t_k))_{L_2(\Gamma)} =(\psi^i_h,u_d(\cdot,T-t_k))_{L_2(\Gamma)}\)
  • \(\boldsymbol{A}\) be the matrices with entries \(\boldsymbol{A}_{ik} = a(s_i,t_k)\)
  • \(\boldsymbol{B}\) be the matrices with entries \(\boldsymbol{B}_{ik} = b(s_i,t_k)\)
  • \(\bar{\boldsymbol{P}} = \text{reversecolumns}(\bar{\boldsymbol{Q}})\)

Observe that the entries of \(\boldsymbol{F}\) are given by \(\boldsymbol{F}_{ik}=(\psi^i_h, f(\cdot,t_k))_{L_2(\Gamma)}\) where \(\psi^i_h\) is the \(i\)-th basis function on a mesh with \(N_h\) nodes. By introducing a finer integration mesh with \(N_{h_{\text{ok}}}\) nodes, we can approximate \(\boldsymbol{F}_{ik}\) as \(\boldsymbol{\psi}_i^\top\boldsymbol{C}^{\text{ok}}\boldsymbol{f}_k\), where \(\boldsymbol{\psi}_i=[\psi^i_h(s_1)\dots \psi^i_h(s_{N_{h_{\text{ok}}}})]^\top\) for \(i=1\dots, N_h\), \(\boldsymbol{f}_k = [f(s_1,t_k)\dots f(s_{N_{h_{\text{ok}}}},t_k)]\top\) for \(k= 0,\dots,N\), and matrix \(\boldsymbol{C}^{\text{ok}}\) has entries \(\boldsymbol{C}^{\text{ok}}_{nm} = (\psi_{h_\text{ok}}^n,\psi_{h_\text{ok}}^m)_{L_2(\Gamma)}\) for \(n,m = 1,\dots, N_{h_\text{ok}}\). Let \(\boldsymbol{f}\) denote the matrix with columns \(\boldsymbol{f}_k\) and \(\boldsymbol{\Psi}\) denote a matrix with entries \(\boldsymbol{\Psi}_{ji}=\psi^i_h(s_j)\) for \(i=1\dots, N_h\) and \(j = 1,\dots N_{h_{\text{ok}}}\). Then \[\begin{equation} \boldsymbol{F} \approx \boldsymbol{\Psi}^\top \boldsymbol{C}^{\text{ok}} \boldsymbol{f}. \end{equation}\] Similarly, if \(\boldsymbol{U_d}\) is the matrix with entries \(\boldsymbol{U_d}^{jk} = u_d(s_j,t_k)\) for \(j = 1,\dots, N_{h_{\text{ok}}}\) and \(k = 0,\dots, N\), then the matrix \(\boldsymbol{V_d}\) can be approximated as \[\begin{equation} \boldsymbol{V_d} \approx \text{reversecolumns}(\boldsymbol{\Psi}^\top \boldsymbol{C}^{\text{ok}} \boldsymbol{U_d}). \end{equation}\]

To approximate \(\boldsymbol{F}\) and \(\boldsymbol{V_d}\), we use the values of \(f\) and \(u_d\) on the integration mesh with \(N_{h_{\text{ok}}}\) nodes. However, to approximate \(\bar{\boldsymbol{Z}}\) and \(\bar{\boldsymbol{V}}\), we do not even have access to the values of \(\bar{z}\) and \(\bar{v}\) on the coarse mesh, let alone on the fine mesh. Instead, we use their approximations on the coarse mesh, given by the matrices \(\bar{\boldsymbol{z}} = \max\{\boldsymbol{A},\min\{\boldsymbol{B},-\bar{\boldsymbol{P}}/\mu\}\}\) and \(\text{reversecolumns}(\bar{\boldsymbol{U}})\), respectively. With these at hand, we can approximate \(\bar{\boldsymbol{Z}}\) and \(\bar{\boldsymbol{V}}\) by first projecting \(\bar{\boldsymbol{z}} = \max\{\boldsymbol{A},\min\{\boldsymbol{B},-\bar{\boldsymbol{P}}/\mu\}\}\) and \(\text{reversecolumns}(\bar{\boldsymbol{U}})\) onto the fine mesh and then approximating the inner products \((\psi^i_h,\bar{z}(\cdot,t_k))_{L_2(\Gamma)}\) and \((\psi^i_h,\bar{u}(\cdot,T-t_k))_{L_2(\Gamma)}\) as \[\begin{equation} \bar{\boldsymbol{Z}} \approx \boldsymbol{\Psi}^\top \boldsymbol{C}^{\text{ok}} \boldsymbol{\Psi} \bar{\boldsymbol{z}}\quad \text{ and }\quad \bar{\boldsymbol{V}} \approx \text{reversecolumns}(\boldsymbol{\Psi}^\top \boldsymbol{C}^{\text{ok}} \boldsymbol{\Psi} \bar{\boldsymbol{U}}). \end{equation}\]

1.3 Construction of an exact solution to the optimal control problem

We follow the approach in Glusa and Otárola (2021) to construct an exact solution to the optimal control problem \(\eqref{eq:costfun}\)-\(\eqref{control_constraints}\). The solution to the elliptic problem \[\begin{align} \label{eq:elliptic_problem} \tag{12} (\kappa^2 - \Delta_\Gamma)^{\alpha/2} u(s)=f(s),\quad s\in\Gamma, \end{align}\] where \(u(\cdot)\) satisfies the Kirchhoff vertex conditions \(\eqref{eq:Kcond}\), is given by \[\begin{align*} u(s) = \sum_{j=1}^{\infty}\lambda_j^{-\alpha/2}(f,e_j)_{L_2(\Gamma)}e_j(s),\quad s\in\Gamma. \end{align*}\] We choose coefficients \(\{x_j\}_{j=1}^{\infty}\) and \(\{y_j\}_{j=1}^{\infty}\) and set \[\begin{align} \label{eq:rhsfunctions} \tag{13} f(s) = \sum_{j=1}^{\infty}x_je_j(s)\quad \text{ and } \quad g(s) = \sum_{j=1}^{\infty}y_je_j(s). \end{align}\] Hence the corresponding solutions are given by \[\begin{align} \label{eq:sols} \tag{14} u(s) = \sum_{j=1}^{\infty}\lambda_j^{-\alpha/2}x_je_j(s)\quad \text{ and } \quad v(s) = \sum_{j=1}^{\infty}\lambda_j^{-\alpha/2}y_je_j(s). \end{align}\] We also define \(\psi(t) = \cos(t)\) and \(\phi(t) =\sin(T-t)\). Observe that \(\psi(0)=1\) and \(\phi(T)=0\). We also choose \(\mu =0.1\), \(a(s,t) = -0.5\) and \(b(s,t)=0.5\). Now we set \[\begin{align*} f(s,t) &= \psi'(t)u(s)+\psi(t)f(s) - \bar{z}(s,t) = -\sin(t)\sum_{j=1}^{\infty}\lambda_j^{-\alpha/2}x_je_j(s)+ \cos(t)\sum_{j=1}^{\infty}x_je_j(s)- \max\{-0.5,\min\{0.5, \sin(T-t)\sum_{j=1}^{\infty}\lambda_j^{-\alpha/2}y_je_j(s)\}\}\\ u_d(s,t) & = \psi(t)u(s)-\mu\phi'(t)v(s)+\mu\phi(t)g(s) = \cos(t)\sum_{j=1}^{\infty}\lambda_j^{-\alpha/2}x_je_j(s)+\mu\cos(T-t)\sum_{j=1}^{\infty}\lambda_j^{-\alpha/2}y_je_j(s)+\mu\sin(T-t)\sum_{j=1}^{\infty}y_je_j(s)\\ u_0(s) & = u(s)= \sum_{j=1}^{\infty}\lambda_j^{-\alpha/2}x_je_j(s) \end{align*}\] The exact solution to the optimal control problem is given by \[\begin{align} \bar{u}(s,t)&=\psi(t)u(s) = \cos(t)\sum_{j=1}^{\infty}\lambda_j^{-\alpha/2}x_je_j(s)\\ \bar{p}(s,t)&=-\mu\phi(t)v(s)= - \mu\sin(T-t)\sum_{j=1}^{\infty}\lambda_j^{-\alpha/2}y_je_j(s)\\ \bar{z}(s,t)&=\max\{a(s,t),\min\{b(s,t), \phi(t)v(s)\}\}= \max\{-0.5,\min\{0.5, \sin(T-t)\sum_{j=1}^{\infty}\lambda_j^{-\alpha/2}y_je_j(s)\}\} \end{align}\]

Observe that the choice of \(f\) decouples the systems \(\eqref{numericalscheme1}\)-\(\eqref{thenumericalscheme2}\).

1.3.1 Checking that \(\bar{u}\) and \(\bar{p}\) satisfy \(\eqref{eq:maineqoptimal}\) and \(\eqref{eq:adjointeq}\)

At \(t=0\), we get \[\begin{align*} \bar{u}(s,0) = \psi(0)u(s) = u(s) = u_0(s), \end{align*}\] while for \((s,t) \in \Gamma \times (0, T)\), we obtain \[\begin{align*} \partial_t \bar{u}(s,t) +(\kappa^2 - \Delta_\Gamma)^{\alpha/2} \bar{u}(s,t) = \psi'(t)u(s) + \psi(t) (\kappa^2 - \Delta_\Gamma)^{\alpha/2} u(s) = \psi'(t)u(s) + \psi(t)f(s) = f(s,t) + \bar{z}(s,t). \end{align*}\] Similarly, at \(t=T\), we get \[\begin{align*} \bar{p}(s,T) = -\mu\phi(T)v(s) = 0, \end{align*}\] while for \((s,t) \in \Gamma \times (0, T)\), we obtain \[\begin{align*} -\partial_t \bar{p}(s,t) +(\kappa^2 - \Delta_\Gamma)^{\alpha/2} \bar{p}(s,t) = \mu\phi'(t)v(s) -\mu \phi(t) (\kappa^2 - \Delta_\Gamma)^{\alpha/2} v(s) = \mu\phi'(t)v(s) -\mu \phi(t)g(s) \end{align*}\] and \[\begin{align*} \bar{u}(s,t) - u_d(s,t) = \psi(t)u(s) - \psi(t)u(s) + \mu\phi'(t)v(s) - \mu\phi(t)g(s) = \mu\phi'(t)v(s) -\mu \phi(t)g(s). \end{align*}\] This shows that \(\bar{u}\) and \(\bar{p}\) satisfy \(\eqref{eq:maineqoptimal}\) and \(\eqref{eq:adjointeq}\).

1.4 Numerical implementation

In this section, we implement the numerical scheme for the optimal control problem \(\eqref{eq:costfun}\)-\(\eqref{control_constraints}\) on the tadpole graph. We start by defining the parameters and constructing the exact solution. We them use \(\eqref{numericalscheme1}\)-\(\eqref{thenumericalscheme2}\) to solve the forward and adjoint problems. Finally, we compute the control variable \(\bar{z}\) using \(\eqref{zz}\).

# Parameters
T_final <- 2
time_step <- 0.02
h <- 0.02
kappa <- 15
alpha <- 1.8
m = 4
beta <- alpha/2
mu <- 0.1
a <- - 0.5
b <- 0.5
N_finite = 4 # choose even
adjusted_N_finite <- N_finite + N_finite/2 + 1
# Coefficients for f and g
coeff_elliptic_g <- rep(0, adjusted_N_finite)
coeff_elliptic_g[7] <- 100
coeff_elliptic_f <- rep(0, adjusted_N_finite)
coeff_elliptic_f[7] <- -10
# time and spatial mesh 
time_seq <- seq(0, T_final, length.out = ((T_final - 0) / time_step + 1))
graph <- gets.graph.tadpole(h = h)
graph$compute_fem()
G <- graph$mesh$G
C <- graph$mesh$C
L <- kappa^2*C + G
I <- Matrix::Diagonal(nrow(C))
# Compute the eigenvalues and eigenfunctions
eigen_params <- gets.eigen.params(N_finite = N_finite, 
                                  kappa = kappa, 
                                  alpha = alpha, 
                                  graph = graph)
EIGENVAL_MINUS_ALPHA <- eigen_params$EIGENVAL_MINUS_ALPHA
EIGENFUN <- eigen_params$EIGENFUN
# Construct the right hand side functions f and g for the elliptic problem
elliptic_f <- as.vector(EIGENFUN %*% coeff_elliptic_f)
elliptic_g <- as.vector(EIGENFUN %*% coeff_elliptic_g)
# Construct the corresponding elliptic solution u and v
elliptic_u <- as.vector(EIGENFUN %*% (coeff_elliptic_f * EIGENVAL_MINUS_ALPHA))
elliptic_v <- as.vector(EIGENFUN %*% (coeff_elliptic_g * EIGENVAL_MINUS_ALPHA))

A <- matrix(a, nrow = nrow(C), ncol = length(time_seq))
B <- matrix(b, nrow = nrow(C), ncol = length(time_seq))
# Construct the optimal variables
psi <- cos(time_seq)
psi_prime <- - sin(time_seq)
phi <- sin(T_final - time_seq)
phi_prime <- - cos(T_final - time_seq)
# psi_rate <- 1
# phi_rate <- 1
# phi_coefficient <- 1
# psi <- exp(psi_rate*time_seq)
# psi_prime <- psi_rate * exp(psi_rate*time_seq)
# phi <- phi_coefficient * exp(phi_rate*time_seq) - phi_coefficient * exp(phi_rate*T_final)
# phi_prime <- phi_coefficient * phi_rate * exp(phi_rate*time_seq)
u_bar <- elliptic_u %*% t(psi)
p_bar <- - mu * (elliptic_v %*% t(phi))
z_bar <- pmax(A, pmin(B, - p_bar / mu))
# Construct the data for the problem
f <- elliptic_u %*% t(psi_prime) + elliptic_f %*% t(psi) - z_bar
f_plus_z_bar <- elliptic_u %*% t(psi_prime) + elliptic_f %*% t(psi)
u_d <- elliptic_u %*% t(psi) -
  mu * (elliptic_v %*% t(phi_prime)) +
  mu * (elliptic_g %*% t(phi))
u_0 <- elliptic_u
# Construct the fractional operator, which is shared for the forward and adjoint problems
my_op_frac <- my.fractional.operators.frac(L, beta, C, scale.factor = kappa^2, m = m, time_step)
# Construct the integration mesh
overkill_graph <- gets.graph.tadpole(h = 0.0001)
overkill_graph$compute_fem()
overkill_EIGENFUN <- gets.eigen.params(N_finite = N_finite, 
                                       kappa = kappa, 
                                       alpha = alpha, 
                                       graph = overkill_graph)$EIGENFUN
# Construct the right hand side functions f and g for the elliptic problem on the integration mesh
overkill_elliptic_f <- as.vector(overkill_EIGENFUN %*% coeff_elliptic_f)
overkill_elliptic_g <- as.vector(overkill_EIGENFUN %*% coeff_elliptic_g)
# Construct the corresponding elliptic solution u and v on the integration mesh
overkill_elliptic_u <- as.vector(overkill_EIGENFUN %*% (coeff_elliptic_f * EIGENVAL_MINUS_ALPHA))
overkill_elliptic_v <- as.vector(overkill_EIGENFUN %*% (coeff_elliptic_g * EIGENVAL_MINUS_ALPHA))
# Construct the projection matrix
Psi <- graph$fem_basis(overkill_graph$get_mesh_locations())
R <- t(Psi) %*% overkill_graph$mesh$C
# Construct the right hand side term for the forward problem
F_plus_Z_bar <- overkill_elliptic_u %*% t(psi_prime) + overkill_elliptic_f %*% t(psi) 
F_plus_Z_bar_innerproduct <- R %*% F_plus_Z_bar
# Solve the forward problem
U_bar <- solve_fractional_evolution(my_op_frac, 
                                     time_step, 
                                     time_seq, 
                                     val_at_0 = u_0, 
                                     RHST = F_plus_Z_bar_innerproduct)
# Construct the right hand side term for the adjoint problem
U_d <- overkill_elliptic_u %*% t(psi) -
  mu * (overkill_elliptic_v %*% t(phi_prime)) + 
  mu * (overkill_elliptic_g %*% t(phi))
V_d <- reversecolumns(R %*% U_d)
V_bar <- reversecolumns(R %*% Psi %*% U_bar)
# Solve the adjoint problem
Q_bar <- solve_fractional_evolution(my_op_frac, 
                                     time_step, 
                                     time_seq, 
                                     val_at_0 = u_0*0, 
                                     RHST = V_bar - V_d)
P_bar <- reversecolumns(Q_bar)
# Compute the control variable
Z_bar <- pmax(A, pmin(B, - P_bar / mu))
start_idx <- which.min(abs(time_seq - 0))
end_idx <- which.min(abs(time_seq - T_final))
idx <- seq(start_idx, end_idx, by = 1)
graph.plotter.3d.static(graph, elliptic_f, elliptic_g)

Figure 1: Right-hand side functions \(f\) and \(g\) for the elliptic problem.

graph.plotter.3d.static(graph, elliptic_u, elliptic_v)

Figure 2: Solutions \(u\) and \(v\) to the elliptic problem.

graph.plotter.3d.single(graph, f_plus_z_bar[, idx], time_seq[idx])

Figure 3: Time evolution of the right-hand side function \(f+\bar{z}\) of the state equation.

v_bar_minus_v_d <- reversecolumns(u_bar - u_d)
graph.plotter.3d.single(graph, v_bar_minus_v_d[, idx], time_seq[idx])

Figure 4: Time evolution of the right-hand side function \(\bar{v}-v_d\) of the adjoint equation.

1.5 Visualization of the results

graph.plotter.3d(graph, time_seq[idx], time_seq[idx], u_bar[, idx], U_bar[, idx])

Figure 5: Exact (in red) and numerical (in blue) optimal state variable.

graph.plotter.3d(graph, time_seq[idx], time_seq[idx], p_bar[, idx], P_bar[, idx])

Figure 6: Exact (in red) and numerical (in blue) optimal adjoint variable.

graph.plotter.3d(graph, time_seq[idx], time_seq[idx], z_bar[, idx], Z_bar[, idx])

Figure 7: Exact (in red) and numerical (in blue) optimal control variable.

graph.plotter.3d(graph, time_seq[idx], time_seq[idx], u_d[, idx], u_bar[, idx])

Figure 8: Desired state (in red) and optimal state (in blue).

2 References

cite_packages(output = "paragraph", out.dir = ".")

We used R version 4.5.0 (R Core Team 2025) and the following R packages: gsignal v. 0.3.7 (Van Boxtel, G.J.M., et al. 2021), here v. 1.0.1 (Müller 2020), htmltools v. 0.5.8.1 (Cheng et al. 2024), knitr v. 1.50 (Xie 2014, 2015, 2025), Matrix v. 1.7.3 (Bates, Maechler, and Jagan 2025), MetricGraph v. 1.5.0.9000 (Bolin, Simas, and Wallin 2023a, 2023b, 2024, 2025; Bolin et al. 2024), patchwork v. 1.3.1 (Pedersen 2025), plotly v. 4.10.4 (Sievert 2020), RColorBrewer v. 1.1.3 (Neuwirth 2022), renv v. 1.0.7 (Ushey and Wickham 2024), reshape2 v. 1.4.4 (Wickham 2007), rmarkdown v. 2.29 (Xie, Allaire, and Grolemund 2018; Xie, Dervieux, and Riederer 2020; Allaire et al. 2024), rSPDE v. 2.5.1.9000 (Bolin and Kirchner 2020; Bolin and Simas 2023; Bolin, Simas, and Xiong 2024), scales v. 1.4.0 (Wickham, Pedersen, and Seidel 2025), slackr v. 3.4.0 (Kaye et al. 2025), tidyverse v. 2.0.0 (Wickham et al. 2019), viridisLite v. 0.4.2 (Garnier et al. 2023), xaringanExtra v. 0.8.0 (Aden-Buie and Warkentin 2024).

Aden-Buie, Garrick, and Matthew T. Warkentin. 2024. xaringanExtra: Extras and Extensions for xaringan Slides. https://doi.org/10.32614/CRAN.package.xaringanExtra.
Allaire, JJ, Yihui Xie, Christophe Dervieux, Jonathan McPherson, Javier Luraschi, Kevin Ushey, Aron Atkins, et al. 2024. rmarkdown: Dynamic Documents for r. https://github.com/rstudio/rmarkdown.
Bates, Douglas, Martin Maechler, and Mikael Jagan. 2025. Matrix: Sparse and Dense Matrix Classes and Methods. https://doi.org/10.32614/CRAN.package.Matrix.
Bolin, David, and Kristin Kirchner. 2020. “The Rational SPDE Approach for Gaussian Random Fields with General Smoothness.” Journal of Computational and Graphical Statistics 29 (2): 274–85. https://doi.org/10.1080/10618600.2019.1665537.
Bolin, David, Mihály Kovács, Vivek Kumar, and Alexandre B. Simas. 2024. “Regularity and Numerical Approximation of Fractional Elliptic Differential Equations on Compact Metric Graphs.” Mathematics of Computation 93 (349): 2439–72. https://doi.org/10.1090/mcom/3929.
Bolin, David, and Alexandre B. Simas. 2023. rSPDE: Rational Approximations of Fractional Stochastic Partial Differential Equations. https://CRAN.R-project.org/package=rSPDE.
Bolin, David, Alexandre B. Simas, and Jonas Wallin. 2023a. MetricGraph: Random Fields on Metric Graphs. https://CRAN.R-project.org/package=MetricGraph.
———. 2023b. “Statistical Inference for Gaussian Whittle-Matérn Fields on Metric Graphs.” arXiv Preprint arXiv:2304.10372. https://doi.org/10.48550/arXiv.2304.10372.
———. 2024. “Gaussian Whittle-Matérn Fields on Metric Graphs.” Bernoulli 30 (2): 1611–39. https://doi.org/10.3150/23-BEJ1647.
———. 2025. “Markov Properties of Gaussian Random Fields on Compact Metric Graphs.” Bernoulli. https://doi.org/10.48550/arXiv.2304.03190.
Bolin, David, Alexandre B. Simas, and Zhen Xiong. 2024. “Covariance-Based Rational Approximations of Fractional SPDEs for Computationally Efficient Bayesian Inference.” Journal of Computational and Graphical Statistics 33 (1): 64–74. https://doi.org/10.1080/10618600.2023.2231051.
Cheng, Joe, Carson Sievert, Barret Schloerke, Winston Chang, Yihui Xie, and Jeff Allen. 2024. htmltools: Tools for HTML. https://doi.org/10.32614/CRAN.package.htmltools.
Garnier, Simon, Ross, Noam, Rudis, Robert, Camargo, et al. 2023. viridis(Lite) - Colorblind-Friendly Color Maps for r. https://doi.org/10.5281/zenodo.4678327.
Glusa, Christian, and Enrique Otárola. 2021. “Error Estimates for the Optimal Control of a Parabolic Fractional PDE.” SIAM J. Numer. Anal. 59 (2): 1140–65. https://doi.org/10.1137/19M1267581.
Kaye, Matt, Bob Rudis, Andrie de Vries, and Jonathan Sidi. 2025. slackr: Send Messages, Images, r Objects and Files to Slack Channels/Users. https://github.com/mrkaye97/slackr.
Müller, Kirill. 2020. here: A Simpler Way to Find Your Files. https://doi.org/10.32614/CRAN.package.here.
Neuwirth, Erich. 2022. RColorBrewer: ColorBrewer Palettes. https://doi.org/10.32614/CRAN.package.RColorBrewer.
Pedersen, Thomas Lin. 2025. patchwork: The Composer of Plots. https://doi.org/10.32614/CRAN.package.patchwork.
R Core Team. 2025. R: A Language and Environment for Statistical Computing. Vienna, Austria: R Foundation for Statistical Computing. https://www.R-project.org/.
Sievert, Carson. 2020. Interactive Web-Based Data Visualization with r, Plotly, and Shiny. Chapman; Hall/CRC. https://plotly-r.com.
Ushey, Kevin, and Hadley Wickham. 2024. renv: Project Environments. https://doi.org/10.32614/CRAN.package.renv.
Van Boxtel, G.J.M., et al. 2021. gsignal: Signal Processing. https://github.com/gjmvanboxtel/gsignal.
Wickham, Hadley. 2007. “Reshaping Data with the reshape Package.” Journal of Statistical Software 21 (12): 1–20. http://www.jstatsoft.org/v21/i12/.
Wickham, Hadley, Mara Averick, Jennifer Bryan, Winston Chang, Lucy D’Agostino McGowan, Romain François, Garrett Grolemund, et al. 2019. “Welcome to the tidyverse.” Journal of Open Source Software 4 (43): 1686. https://doi.org/10.21105/joss.01686.
Wickham, Hadley, Thomas Lin Pedersen, and Dana Seidel. 2025. scales: Scale Functions for Visualization. https://doi.org/10.32614/CRAN.package.scales.
Xie, Yihui. 2014. knitr: A Comprehensive Tool for Reproducible Research in R.” In Implementing Reproducible Computational Research, edited by Victoria Stodden, Friedrich Leisch, and Roger D. Peng. Chapman; Hall/CRC.
———. 2015. Dynamic Documents with R and Knitr. 2nd ed. Boca Raton, Florida: Chapman; Hall/CRC. https://yihui.org/knitr/.
———. 2025. knitr: A General-Purpose Package for Dynamic Report Generation in R. https://yihui.org/knitr/.
Xie, Yihui, J. J. Allaire, and Garrett Grolemund. 2018. R Markdown: The Definitive Guide. Boca Raton, Florida: Chapman; Hall/CRC. https://bookdown.org/yihui/rmarkdown.
Xie, Yihui, Christophe Dervieux, and Emily Riederer. 2020. R Markdown Cookbook. Boca Raton, Florida: Chapman; Hall/CRC. https://bookdown.org/yihui/rmarkdown-cookbook.
LS0tCnRpdGxlOiAiT3B0aW1hbCBjb250cm9sIG9mIGZyYWN0aW9uYWwgZGlmZnVzaW9uIGVxdWF0aW9ucyBvbiBtZXRyaWMgZ3JhcGhzIgpkYXRlOiAiTGFzdCBtb2RpZmllZDogYHIgZm9ybWF0KFN5cy50aW1lKCksICclZC0lbS0lWS4nKWAiCm91dHB1dDoKICBodG1sX2RvY3VtZW50OgogICAgbWF0aGpheDogImh0dHBzOi8vY2RuLmpzZGVsaXZyLm5ldC9ucG0vbWF0aGpheEAzL2VzNS90ZXgtbW1sLWNodG1sLmpzIgogICAgaGlnaGxpZ2h0OiBweWdtZW50cwogICAgdGhlbWU6IGZsYXRseQogICAgY29kZV9mb2xkaW5nOiBoaWRlICMgY2xhc3Muc291cmNlID0gImZvbGQtaGlkZSIgdG8gaGlkZSBjb2RlIGFuZCBhZGQgYSBidXR0b24gdG8gc2hvdyBpdAogICAgZGZfcHJpbnQ6IHBhZ2VkCiAgICB0b2M6IHRydWUKICAgIHRvY19mbG9hdDoKICAgICAgY29sbGFwc2VkOiB0cnVlCiAgICAgIHNtb290aF9zY3JvbGw6IHRydWUKICAgIG51bWJlcl9zZWN0aW9uczogdHJ1ZQogICAgZmlnX2NhcHRpb246IHRydWUKICAgIGNvZGVfZG93bmxvYWQ6IHRydWUKICAgIGNzczogdmlzdWFsLmNzcwphbHdheXNfYWxsb3dfaHRtbDogdHJ1ZQpiaWJsaW9ncmFwaHk6IAogIC0gcmVmZXJlbmNlcy5iaWIKICAtIGdyYXRlZnVsLXJlZnMuYmliCmhlYWRlci1pbmNsdWRlczoKICAtIFxuZXdjb21tYW5ke1xhcn17XG1hdGhiYntSfX0KICAtIFxuZXdjb21tYW5ke1xsbGF2fVsxXXtcbGVmdFx7IzFccmlnaHRcfX0KICAtIFxuZXdjb21tYW5ke1xwYXJlfVsxXXtcbGVmdCgjMVxyaWdodCl9CiAgLSBcbmV3Y29tbWFuZHtcTmNhbH17XG1hdGhjYWx7Tn19CiAgLSBcbmV3Y29tbWFuZHtcVmNhbH17XG1hdGhjYWx7Vn19CiAgLSBcbmV3Y29tbWFuZHtcRWNhbH17XG1hdGhjYWx7RX19CiAgLSBcbmV3Y29tbWFuZHtcV2NhbH17XG1hdGhjYWx7V319CiAgLSBcbmV3Y29tbWFuZHtcYWxtb3N0ZXZlcnl3aGVyZX17XG1hdGhybXthLmUufVw7fQotLS0KCkdvIGJhY2sgdG8gdGhlIFtDb250ZW50c10oYWJvdXQuaHRtbCkgcGFnZS4KCjxkaXYgc3R5bGU9ImNvbG9yOiAjMmMzZTUwOyB0ZXh0LWFsaWduOiByaWdodDsiPgoqKioqKioqKiAgCjxzdHJvbmc+UHJlc3MgU2hvdyB0byByZXZlYWwgdGhlIGNvZGUgY2h1bmtzLjwvc3Ryb25nPiAgCgoqKioqKioqKgo8L2Rpdj4KCgpgYGB7cn0KIyBDcmVhdGUgYSBjbGlwYm9hcmQgYnV0dG9uIG9uIHRoZSByZW5kZXJlZCBIVE1MIHBhZ2UKc291cmNlKGhlcmU6OmhlcmUoImNsaXBib2FyZC5SIikpOyBjbGlwYm9hcmQKIyBTZXQgc2VlZCBmb3IgcmVwcm9kdWNpYmlsaXR5CnNldC5zZWVkKDE5ODIpIAojIFNldCBnbG9iYWwgb3B0aW9ucyBmb3IgYWxsIGNvZGUgY2h1bmtzCmtuaXRyOjpvcHRzX2NodW5rJHNldCgKICAjIERpc2FibGUgbWVzc2FnZXMgcHJpbnRlZCBieSBSIGNvZGUgY2h1bmtzCiAgbWVzc2FnZSA9IEZBTFNFLCAgICAKICAjIERpc2FibGUgd2FybmluZ3MgcHJpbnRlZCBieSBSIGNvZGUgY2h1bmtzCiAgd2FybmluZyA9IEZBTFNFLCAgICAKICAjIFNob3cgUiBjb2RlIHdpdGhpbiBjb2RlIGNodW5rcyBpbiBvdXRwdXQKICBlY2hvID0gVFJVRSwgICAgICAgIAogICMgSW5jbHVkZSBib3RoIFIgY29kZSBhbmQgaXRzIHJlc3VsdHMgaW4gb3V0cHV0CiAgaW5jbHVkZSA9IFRSVUUsICAgICAKICAjIEV2YWx1YXRlIFIgY29kZSBjaHVua3MKICBldmFsID0gVFJVRSwgICAgICAgCiAgIyBFbmFibGUgY2FjaGluZyBvZiBSIGNvZGUgY2h1bmtzIGZvciBmYXN0ZXIgcmVuZGVyaW5nCiAgY2FjaGUgPSBGQUxTRSwgICAgICAKICAjIEFsaWduIGZpZ3VyZXMgaW4gdGhlIGNlbnRlciBvZiB0aGUgb3V0cHV0CiAgZmlnLmFsaWduID0gImNlbnRlciIsCiAgIyBFbmFibGUgcmV0aW5hIGRpc3BsYXkgZm9yIGhpZ2gtcmVzb2x1dGlvbiBmaWd1cmVzCiAgcmV0aW5hID0gMiwKICAjIFNob3cgZXJyb3JzIGluIHRoZSBvdXRwdXQgaW5zdGVhZCBvZiBzdG9wcGluZyByZW5kZXJpbmcKICBlcnJvciA9IFRSVUUsCiAgIyBEbyBub3QgY29sbGFwc2UgY29kZSBhbmQgb3V0cHV0IGludG8gYSBzaW5nbGUgYmxvY2sKICBjb2xsYXBzZSA9IEZBTFNFCikKIyBTdGFydCB0aGUgZmlndXJlIGNvdW50ZXIKZmlnX2NvdW50IDwtIDAKIyBEZWZpbmUgdGhlIGNhcHRpb25lciBmdW5jdGlvbgpjYXB0aW9uZXIgPC0gZnVuY3Rpb24oY2FwdGlvbikgewogIGZpZ19jb3VudCA8PC0gZmlnX2NvdW50ICsgMQogIHBhc3RlMCgiRmlndXJlICIsIGZpZ19jb3VudCwgIjogIiwgY2FwdGlvbikKfQpgYGAKCmBgYHtyfQpsaWJyYXJ5KE1ldHJpY0dyYXBoKQpsaWJyYXJ5KGdncGxvdDIpCmxpYnJhcnkocmVzaGFwZTIpCmxpYnJhcnkocGxvdGx5KQpgYGAKCgpgYGB7cn0KY2FwdHVyZS5vdXRwdXQoCiAga25pdHI6OnB1cmwoaGVyZTo6aGVyZSgiZnVuY3Rpb25hbGl0eS5SbWQiKSwgb3V0cHV0ID0gaGVyZTo6aGVyZSgiZnVuY3Rpb25hbGl0eS5SIikpLAogIGZpbGUgPSBoZXJlOjpoZXJlKCJwdXJsX2xvZy50eHQiKQopCnNvdXJjZShoZXJlOjpoZXJlKCJmdW5jdGlvbmFsaXR5LlIiKSkKYGBgCgojIE9wdGltYWwgY29udHJvbCBvZiBmcmFjdGlvbmFsIGRpZmZ1c2lvbiBlcXVhdGlvbnMgb24gbWV0cmljIGdyYXBocwoKIyMgUHJvYmxlbSBzdGF0ZW1lbnQgeyNvcHRpbWFsX2NvbnRyb2xfcHJvYmxlbX0KCkxldCAkXEdhbW1hID0gKFxtYXRoY2Fse1Z9LFxtYXRoY2Fse0V9KSQgYmUgYSBtZXRyaWMgZ3JhcGguIExldCAkdV9kOiBcR2FtbWEgXHRpbWVzKDAsIFQpIFxyaWdodGFycm93IFxtYXRoYmJ7Un0kIGJlIHRoZSBkZXNpcmVkIHN0YXRlIGFuZCAgJFxtdT4wJCBhIHJlZ3VsYXJpemF0aW9uIHBhcmFtZXRlci4gV2UgZGVmaW5lIHRoZSBjb3N0IGZ1bmN0aW9uYWwKClxiZWdpbntlcXVhdGlvbn0KXGxhYmVse2VxOmNvc3RmdW59Clx0YWd7MX0KICAgIEoodSwgeik9XGZyYWN7MX17Mn0gXGludF8wXlRcbGVmdChcbGVmdFx8dS11X2RccmlnaHRcfF97TF8yKFxHYW1tYSl9XjIrXG11XHx6XHxfe0xfMihcR2FtbWEpfV4yXHJpZ2h0KSBkdApcZW5ke2VxdWF0aW9ufQoKCkxldCAkZjpcR2FtbWFcdGltZXMgKDAsVClccmlnaHRhcnJvd1xhciQgYW5kICR1XzA6IFxHYW1tYSBccmlnaHRhcnJvdyBcbWF0aGJie1J9JCBiZSBmaXhlZCBmdW5jdGlvbnMuIFdlIHdpbGwgY2FsbCB0aGVtIHJpZ2h0LWhhbmQgc2lkZSBhbmQgaW5pdGlhbCBkYXR1bSwgcmVzcGVjdGl2ZWx5LiBMZXQgJFxhbHBoYVxpbigwLDJdJCBhbmQgJHo6IFxHYW1tYSBcdGltZXMoMCwgVCkgXHJpZ2h0YXJyb3cgXG1hdGhiYntSfSQgZGVub3RlIHRoZSBjb250cm9sIHZhcmlhYmxlLiBXZSBzaGFsbCBiZSBjb25jZXJuZWQgd2l0aCB0aGUgZm9sbG93aW5nIFBERS1jb25zdHJhaW5lZCBvcHRpbWl6YXRpb24gcHJvYmxlbTogRmluZAoKXGJlZ2lue2VxdWF0aW9ufQpcbGFiZWx7ZXE6bWluX3Byb30KXHRhZ3syfQogICAgXG1pblw7IEoodSwgeikKXGVuZHtlcXVhdGlvbn0Kc3ViamVjdCB0byB0aGUgZnJhY3Rpb25hbCBkaWZmdXNpb24gZXF1YXRpb24KXGJlZ2lue2VxdWF0aW9ufQpcbGFiZWx7ZXE6bWFpbmVxfQpcdGFnezN9ClxsZWZ0XHsKXGJlZ2lue2FsaWduZWR9CiAgICBccGFydGlhbF90IHUocyx0KSArIChca2FwcGFeMiAtIFxEZWx0YV9cR2FtbWEpXntcYWxwaGEvMn0gdShzLHQpICY9IGYocyx0KSt6KHMsdCksICYmIFxxdWFkIChzLHQpIFxpbiBcR2FtbWEgXHRpbWVzICgwLCBUKSwgXFwKICAgIHUocywwKSAmPSB1XzAocyksICYmIFxxdWFkIHMgXGluIFxHYW1tYSwKXGVuZHthbGlnbmVkfQpccmlnaHQuClxlbmR7ZXF1YXRpb259CndpdGggJHUoXGNkb3QsdCkkIHNhdGlzZnlpbmcgdGhlIEtpcmNoaG9mZiB2ZXJ0ZXggY29uZGl0aW9ucwpcYmVnaW57ZXF1YXRpb259ClxsYWJlbHtlcTpLY29uZH0KXHRhZ3s0fQogICBcbWF0aGNhbHtLfSA9ICBcbGVmdFx7XHBoaVxpbiBDKFxHYW1tYSlcO1xtaWRkbGV8XDsgXGZvcmFsbCB2XGluIFxtYXRoY2Fse1Z9Olw7IFxzdW1fe2VcaW5cbWF0aGNhbHtFfV92fVxwYXJ0aWFsX2UgXHBoaSh2KT0wIFxyaWdodFx9ClxlbmR7ZXF1YXRpb259CmFuZCB0aGUgY29udHJvbCBjb25zdHJhaW50cyAKXGJlZ2lue2FsaWdufQpcbGFiZWx7Y29udHJvbF9jb25zdHJhaW50c30KXHRhZ3s1fQogICAgYShzLHQpXGxlcSB6KHMsdClcbGVxIGIocyx0KVw7XGFsbW9zdGV2ZXJ5d2hlcmUgKHMsdClcaW5cR2FtbWEgXHRpbWVzKDAsIFQpLgpcZW5ke2FsaWdufQoKCiMjIE51bWVyaWNhbCBzY2hlbWUgeyNudW1fc2NoZW1lX29wdGltX2NvbnRyb2x9CgpUaGUgb3B0aW1hbCB2YXJpYWJsZXMgJChcYmFye3V9LCBcYmFye3B9LCBcYmFye3p9KSQgc2F0aXNmeQoKXGJlZ2lue2VxdWF0aW9ufQpcbGFiZWx7ZXE6bWFpbmVxb3B0aW1hbH0KXHRhZ3s2fQpcbGVmdFx7ClxiZWdpbnthbGlnbmVkfQogICAgXHBhcnRpYWxfdCBcYmFye3V9KHMsdCkgKyAoXGthcHBhXjIgLSBcRGVsdGFfXEdhbW1hKV57XGFscGhhLzJ9IFxiYXJ7dX0ocyx0KSAmPSBmKHMsdCkrXGJhcnt6fShzLHQpLCAmJiBccXVhZCAocyx0KSBcaW4gXEdhbW1hIFx0aW1lcyAoMCwgVCksIFxcCiAgICBcYmFye3V9KHMsMCkgJj0gdV8wKHMpLCAmJiBccXVhZCBzIFxpbiBcR2FtbWEsClxlbmR7YWxpZ25lZH0KXHJpZ2h0LgpcZW5ke2VxdWF0aW9ufQphbmQKXGJlZ2lue2VxdWF0aW9ufQpcbGFiZWx7ZXE6YWRqb2ludGVxfQpcdGFnezd9ClxsZWZ0XHsKXGJlZ2lue2FsaWduZWR9CiAgICAtXHBhcnRpYWxfdCBcYmFye3B9KHMsdCkgKyAoXGthcHBhXjIgLSBcRGVsdGFfXEdhbW1hKV57XGFscGhhLzJ9IFxiYXJ7cH0ocyx0KSAmPSBcYmFye3V9KHMsdCktdV9kKHMsdCksICYmIFxxdWFkIChzLHQpIFxpbiBcR2FtbWEgXHRpbWVzICgwLCBUKSwgXFwKICAgIFxiYXJ7cH0ocyxUKSAmPSAwLCAmJiBccXVhZCBzIFxpbiBcR2FtbWEsClxlbmR7YWxpZ25lZH0KXHJpZ2h0LgpcZW5ke2VxdWF0aW9ufQp3aXRoIApcYmVnaW57YWxpZ259ClxsYWJlbHt6en0KXHRhZ3s4fQogICAgICAgIFxiYXJ7en0ocyx0KSA9IFxtYXhcbGVmdFx7YShzLHQpLFxtaW5cbGVmdFx7YihzLHQpLC1cZGZyYWN7MX17XG11fVxiYXJ7cH0ocyx0KVxyaWdodFx9XHJpZ2h0XH0uClxlbmR7YWxpZ259CkJ5IGNvbnNpZGVyaW5nIHRoZSBjaGFuZ2Ugb2YgdmFyaWFibGUgJFxiYXJ7cX0ocyx0KSA9IFxiYXJ7cH0ocyxULXQpJCwgdGhlIGZyYWN0aW9uYWwgYWRqb2ludCBwcm9ibGVtIFxlcXJlZntlcTphZGpvaW50ZXF9IGJlY29tZXMgYSBmb3J3YXJkLWluLXRpbWUgcHJvYmxlbSB3aGVyZSB0aGUgdHJhbnNmb3JtZWQgYWRqb2ludCBzdGF0ZSAkXGJhcntxfSQgc29sdmVzClxiZWdpbntlcXVhdGlvbn0KXGxhYmVse3RyYW5zZm9ybWVkX2Fkam9pbnRfc3RhdGV9Clx0YWd7OX0KXGxlZnRcewpcYmVnaW57YWxpZ25lZH0KICAgIFxwYXJ0aWFsX3QgXGJhcntxfShzLHQpICsgKFxrYXBwYV4yIC0gXERlbHRhX1xHYW1tYSlee1xhbHBoYS8yfSBcYmFye3F9KHMsdCkgJj0gXGJhcnt2fShzLHQpLXZfZChzLHQpLCAmJiBccXVhZCAocyx0KSBcaW4gXEdhbW1hIFx0aW1lcyAoMCwgVCksIFxcCiAgICBcYmFye3F9KHMsMCkgJj0gMCwgJiYgXHF1YWQgcyBcaW4gXEdhbW1hLApcZW5ke2FsaWduZWR9ClxyaWdodC4KXGVuZHtlcXVhdGlvbn0Kc2luY2UgJFxwYXJ0aWFsX3RcYmFye3F9KHMsdCkgPSAtXHBhcnRpYWxfdFxiYXJ7cH0ocyxULXQpJCBhbmQgJFxiYXJ7cX0ocywgMCkgPVxiYXJ7cH0ocyxULTApPSBcYmFye3B9KHMsVCk9MCQuIEhlcmUsICRcYmFye3Z9KHMsdCkgPSBcYmFye3V9KHMsVC10KSQgYW5kICR2X2Qocyx0KSA9IHVfZChzLFQtdCkkLiAKCgpHaXZlbiAkXGJhcnt1fSQgYW5kICR1X2QkLCB3ZSBjYW4gdGltZS1yZXZlcnNlIHRoZW0gdG8gb2J0YWluICRcYmFye3Z9JCBhbmQgJHZfZCQgYW5kIHRoZW4gdXNlIHRoZSBzYW1lIG51bWVyaWNhbCBzY2hlbWUgd2UgdXNlIGZvciB0aGUgZm9yd2FyZCBwcm9ibGVtIFxlcXJlZntlcTptYWluZXFvcHRpbWFsfSB0byBzb2x2ZSB0aGUgYWRqb2ludCBwcm9ibGVtIFxlcXJlZnt0cmFuc2Zvcm1lZF9hZGpvaW50X3N0YXRlfS4gVGhlIGNvbnRyb2wgdmFyaWFibGUgJFxiYXJ7en0kIGlzIHRoZW4gY29tcHV0ZWQgdXNpbmcgXGVxcmVme3p6fS4KClRoZSBudW1lcmljYWwgc2NoZW1lIGZvciBcZXFyZWZ7ZXE6bWFpbmVxb3B0aW1hbH0gYW5kIFxlcXJlZnt0cmFuc2Zvcm1lZF9hZGpvaW50X3N0YXRlfSBhcmUgZ2l2ZW4gYnkgKHNlZSB0aGUgW0Z1bmN0aW9uYWxpdHldKGZ1bmN0aW9uYWxpdHkuaHRtbCNudW1fc2NoZW1lKSBwYWdlKQpcYmVnaW57YWxpZ259ClxsYWJlbHtudW1lcmljYWxzY2hlbWUxfQpcdGFnezEwfQpcYmFye1xib2xkc3ltYm9se1V9fV57aysxfSA9IFxib2xkc3ltYm9se0N9XnstMX1cbGVmdChcc3VtX3trPTF9XnttKzF9IGFfa1xsZWZ0KFxkZnJhY3tcYm9sZHN5bWJvbHtMfVxib2xkc3ltYm9se0N9XnstMX19e1xrYXBwYV4yfS1wX2tcYm9sZHN5bWJvbHtJfVxyaWdodCleey0xfSArIHJcYm9sZHN5bWJvbHtJfVxyaWdodCkgKFxib2xkc3ltYm9se0N9XGJhcntcYm9sZHN5bWJvbHtVfX1eaytcdGF1IChcYm9sZHN5bWJvbHtGfV57aysxfStcYmFye1xib2xkc3ltYm9se1p9fV57aysxfSkpLApcZW5ke2FsaWdufQphbmQKXGJlZ2lue2FsaWdufQpcbGFiZWx7dGhlbnVtZXJpY2Fsc2NoZW1lMn0KXHRhZ3sxMX0KXGJhcntcYm9sZHN5bWJvbHtRfX1ee2srMX0gPSBcYm9sZHN5bWJvbHtDfV57LTF9XGxlZnQoXHN1bV97az0xfV57bSsxfSBhX2tcbGVmdChcZGZyYWN7XGJvbGRzeW1ib2x7TH1cYm9sZHN5bWJvbHtDfV57LTF9fXtca2FwcGFeMn0tcF9rXGJvbGRzeW1ib2x7SX1ccmlnaHQpXnstMX0gKyByXGJvbGRzeW1ib2x7SX1ccmlnaHQpIChcYm9sZHN5bWJvbHtDfVxiYXJ7XGJvbGRzeW1ib2x7UX19XmsrXHRhdSAoXGJhcntcYm9sZHN5bWJvbHtWfX1ee2srMX0tXGJvbGRzeW1ib2x7Vl9kfV57aysxfSkpLApcZW5ke2FsaWdufQp3aXRoIGluaXRpYWwgY29uZGl0aW9ucyAkXGJhcntcYm9sZHN5bWJvbHtVfX1eMCA9IFxib2xkc3ltYm9se1V9XzAkIGFuZCAkXGJhcntcYm9sZHN5bWJvbHtRfX1eMCA9IFxib2xkc3ltYm9sezB9JC4gT2JzZXJ2ZSB0aGF0IFxlcXJlZntudW1lcmljYWxzY2hlbWUxfS1cZXFyZWZ7dGhlbnVtZXJpY2Fsc2NoZW1lMn0gaXMgYSBjb3VwbGVkIHByb2JsZW0uCgpGb3Igbm90YXRpb24gY29udmVuaWVuY2UsIGxldAoKIC0gJFxiYXJ7XGJvbGRzeW1ib2x7VX19JCAgYmUgdGhlIG1hdHJpeCB3aXRoIGNvbHVtbnMgJFxiYXJ7XGJvbGRzeW1ib2x7VX19XmskIG9idGFpbmVkIGZyb20gdGhlIG51bWVyaWNhbCBzY2hlbWUgXGVxcmVme251bWVyaWNhbHNjaGVtZTF9IAogLSAkXGJvbGRzeW1ib2x7Rn0kIGJlIHRoZSBtYXRyaXggd2l0aCBjb2x1bW5zICRcYm9sZHN5bWJvbHtGfV5rJCBhbmQgZW50cmllcyAkXGJvbGRzeW1ib2x7Rn1fe2lrfT0oXHBzaV5pX2gsZihcY2RvdCx0X2spKV97TF8yKFxHYW1tYSl9JAogLSAkXGJhcntcYm9sZHN5bWJvbHtafX0kIGJlIHRoZSBtYXRyaXggd2l0aCBjb2x1bW5zICRcYmFye1xib2xkc3ltYm9se1p9fV5rJCBhbmQgZW50cmllcyAkXGJhcntcYm9sZHN5bWJvbHtafX1fe2lrfT0oXHBzaV5pX2gsXGJhcnt6fShcY2RvdCx0X2spKV97TF8yKFxHYW1tYSl9JAogLSAkXGJhcntcYm9sZHN5bWJvbHtRfX0kIGJlIHRoZSBtYXRyaXggd2l0aCBjb2x1bW5zICRcYmFye1xib2xkc3ltYm9se1F9fV5rJCBvYnRhaW5lZCBmcm9tIHRoZSBudW1lcmljYWwgc2NoZW1lIFxlcXJlZnt0aGVudW1lcmljYWxzY2hlbWUyfQogLSAkXGJhcntcYm9sZHN5bWJvbHtWfX0kIGJlIHRoZSBtYXRyaXggd2l0aCBjb2x1bW5zICRcYmFye1xib2xkc3ltYm9se1Z9fV5rJCBhbmQgZW50cmllcyAkXGJhcntcYm9sZHN5bWJvbHtWfX1ee2lrfT0oXHBzaV5pX2gsXGJhcnt2fShcY2RvdCx0X2spKV97TF8yKFxHYW1tYSl9ID0oXHBzaV5pX2gsXGJhcnt1fShcY2RvdCxULXRfaykpX3tMXzIoXEdhbW1hKX0kCiAtICRcYm9sZHN5bWJvbHtWX2R9JCBiZSB0aGUgbWF0cml4IHdpdGggY29sdW1ucyAkXGJvbGRzeW1ib2x7Vl9kfV5rJCBhbmQgZW50cmllcyAkXGJvbGRzeW1ib2x7Vl9kfV57aWt9PShccHNpXmlfaCx2X2QoXGNkb3QsdF9rKSlfe0xfMihcR2FtbWEpfSA9KFxwc2leaV9oLHVfZChcY2RvdCxULXRfaykpX3tMXzIoXEdhbW1hKX0kCiAtICRcYm9sZHN5bWJvbHtBfSQgYmUgdGhlIG1hdHJpY2VzIHdpdGggZW50cmllcyAkXGJvbGRzeW1ib2x7QX1fe2lrfSA9IGEoc19pLHRfaykkCiAtICRcYm9sZHN5bWJvbHtCfSQgYmUgdGhlIG1hdHJpY2VzIHdpdGggZW50cmllcyAkXGJvbGRzeW1ib2x7Qn1fe2lrfSA9IGIoc19pLHRfaykkCiAtICRcYmFye1xib2xkc3ltYm9se1B9fSA9IFx0ZXh0e3JldmVyc2Vjb2x1bW5zfShcYmFye1xib2xkc3ltYm9se1F9fSkkCgoKT2JzZXJ2ZSB0aGF0IHRoZSBlbnRyaWVzIG9mICRcYm9sZHN5bWJvbHtGfSQgYXJlIGdpdmVuIGJ5ICRcYm9sZHN5bWJvbHtGfV97aWt9PShccHNpXmlfaCwgZihcY2RvdCx0X2spKV97TF8yKFxHYW1tYSl9JCB3aGVyZSAkXHBzaV5pX2gkIGlzIHRoZSAkaSQtdGggYmFzaXMgZnVuY3Rpb24gb24gYSBtZXNoIHdpdGggJE5faCQgbm9kZXMuIEJ5IGludHJvZHVjaW5nIGEgZmluZXIgaW50ZWdyYXRpb24gbWVzaCB3aXRoICROX3toX3tcdGV4dHtva319fSQgbm9kZXMsIHdlIGNhbiBhcHByb3hpbWF0ZSAkXGJvbGRzeW1ib2x7Rn1fe2lrfSQgYXMgJFxib2xkc3ltYm9se1xwc2l9X2leXHRvcFxib2xkc3ltYm9se0N9XntcdGV4dHtva319XGJvbGRzeW1ib2x7Zn1fayQsIHdoZXJlICRcYm9sZHN5bWJvbHtccHNpfV9pPVtccHNpXmlfaChzXzEpXGRvdHMgXHBzaV5pX2goc197Tl97aF97XHRleHR7b2t9fX19KV1eXHRvcCQgZm9yICRpPTFcZG90cywgTl9oJCwgJFxib2xkc3ltYm9se2Z9X2sgPSBbZihzXzEsdF9rKVxkb3RzIGYoc197Tl97aF97XHRleHR7b2t9fX19LHRfayldXHRvcCQgZm9yICRrPSAwLFxkb3RzLE4kLCBhbmQgbWF0cml4ICRcYm9sZHN5bWJvbHtDfV57XHRleHR7b2t9fSQgaGFzIGVudHJpZXMgJFxib2xkc3ltYm9se0N9XntcdGV4dHtva319X3tubX0gPSAoXHBzaV97aF9cdGV4dHtva319Xm4sXHBzaV97aF9cdGV4dHtva319Xm0pX3tMXzIoXEdhbW1hKX0kIGZvciAkbixtID0gMSxcZG90cywgTl97aF9cdGV4dHtva319JC4gTGV0ICRcYm9sZHN5bWJvbHtmfSQgZGVub3RlIHRoZSBtYXRyaXggd2l0aCBjb2x1bW5zICRcYm9sZHN5bWJvbHtmfV9rJCBhbmQgJFxib2xkc3ltYm9se1xQc2l9JCBkZW5vdGUgYSBtYXRyaXggd2l0aCBlbnRyaWVzICRcYm9sZHN5bWJvbHtcUHNpfV97aml9PVxwc2leaV9oKHNfaikkIGZvciAkaT0xXGRvdHMsIE5faCQgYW5kICRqID0gMSxcZG90cyBOX3toX3tcdGV4dHtva319fSQuIFRoZW4gClxiZWdpbntlcXVhdGlvbn0KXGJvbGRzeW1ib2x7Rn0gXGFwcHJveCBcYm9sZHN5bWJvbHtcUHNpfV5cdG9wIFxib2xkc3ltYm9se0N9XntcdGV4dHtva319IFxib2xkc3ltYm9se2Z9LgpcZW5ke2VxdWF0aW9ufQpTaW1pbGFybHksIGlmICRcYm9sZHN5bWJvbHtVX2R9JCBpcyB0aGUgbWF0cml4IHdpdGggZW50cmllcyAkXGJvbGRzeW1ib2x7VV9kfV57amt9ID0gdV9kKHNfaix0X2spJCBmb3IgJGogPSAxLFxkb3RzLCBOX3toX3tcdGV4dHtva319fSQgYW5kICRrID0gMCxcZG90cywgTiQsIHRoZW4gdGhlIG1hdHJpeCAkXGJvbGRzeW1ib2x7Vl9kfSQgY2FuIGJlIGFwcHJveGltYXRlZCBhcwpcYmVnaW57ZXF1YXRpb259Clxib2xkc3ltYm9se1ZfZH0gXGFwcHJveCBcdGV4dHtyZXZlcnNlY29sdW1uc30oXGJvbGRzeW1ib2x7XFBzaX1eXHRvcCBcYm9sZHN5bWJvbHtDfV57XHRleHR7b2t9fSBcYm9sZHN5bWJvbHtVX2R9KS4KXGVuZHtlcXVhdGlvbn0KCgpUbyBhcHByb3hpbWF0ZSAkXGJvbGRzeW1ib2x7Rn0kIGFuZCAkXGJvbGRzeW1ib2x7Vl9kfSQsIHdlIHVzZSB0aGUgdmFsdWVzIG9mICRmJCBhbmQgJHVfZCQgb24gdGhlIGludGVncmF0aW9uIG1lc2ggd2l0aCAkTl97aF97XHRleHR7b2t9fX0kIG5vZGVzLiBIb3dldmVyLCB0byBhcHByb3hpbWF0ZSAkXGJhcntcYm9sZHN5bWJvbHtafX0kIGFuZCAkXGJhcntcYm9sZHN5bWJvbHtWfX0kLCB3ZSBkbyBub3QgZXZlbiBoYXZlIGFjY2VzcyB0byB0aGUgdmFsdWVzIG9mICRcYmFye3p9JCBhbmQgJFxiYXJ7dn0kIG9uIHRoZSBjb2Fyc2UgbWVzaCwgbGV0IGFsb25lIG9uIHRoZSBmaW5lIG1lc2guIEluc3RlYWQsIHdlIHVzZSB0aGVpciBhcHByb3hpbWF0aW9ucyBvbiB0aGUgY29hcnNlIG1lc2gsIGdpdmVuIGJ5IHRoZSBtYXRyaWNlcyAkXGJhcntcYm9sZHN5bWJvbHt6fX0gPSBcbWF4XHtcYm9sZHN5bWJvbHtBfSxcbWluXHtcYm9sZHN5bWJvbHtCfSwtXGJhcntcYm9sZHN5bWJvbHtQfX0vXG11XH1cfSQgYW5kICRcdGV4dHtyZXZlcnNlY29sdW1uc30oXGJhcntcYm9sZHN5bWJvbHtVfX0pJCwgcmVzcGVjdGl2ZWx5LiBXaXRoIHRoZXNlIGF0IGhhbmQsIAp3ZSBjYW4gYXBwcm94aW1hdGUgJFxiYXJ7XGJvbGRzeW1ib2x7Wn19JCBhbmQgJFxiYXJ7XGJvbGRzeW1ib2x7Vn19JCBieSBmaXJzdCBwcm9qZWN0aW5nICRcYmFye1xib2xkc3ltYm9se3p9fSA9IFxtYXhce1xib2xkc3ltYm9se0F9LFxtaW5ce1xib2xkc3ltYm9se0J9LC1cYmFye1xib2xkc3ltYm9se1B9fS9cbXVcfVx9JCBhbmQgJFx0ZXh0e3JldmVyc2Vjb2x1bW5zfShcYmFye1xib2xkc3ltYm9se1V9fSkkIG9udG8gdGhlIGZpbmUgbWVzaCBhbmQgdGhlbiBhcHByb3hpbWF0aW5nIHRoZSBpbm5lciBwcm9kdWN0cyAkKFxwc2leaV9oLFxiYXJ7en0oXGNkb3QsdF9rKSlfe0xfMihcR2FtbWEpfSQgYW5kICQoXHBzaV5pX2gsXGJhcnt1fShcY2RvdCxULXRfaykpX3tMXzIoXEdhbW1hKX0kIGFzClxiZWdpbntlcXVhdGlvbn0KXGJhcntcYm9sZHN5bWJvbHtafX0gXGFwcHJveCBcYm9sZHN5bWJvbHtcUHNpfV5cdG9wIFxib2xkc3ltYm9se0N9XntcdGV4dHtva319IFxib2xkc3ltYm9se1xQc2l9IFxiYXJ7XGJvbGRzeW1ib2x7en19XHF1YWQgXHRleHR7IGFuZCB9XHF1YWQgXGJhcntcYm9sZHN5bWJvbHtWfX0gXGFwcHJveCBcdGV4dHtyZXZlcnNlY29sdW1uc30oXGJvbGRzeW1ib2x7XFBzaX1eXHRvcCBcYm9sZHN5bWJvbHtDfV57XHRleHR7b2t9fSBcYm9sZHN5bWJvbHtcUHNpfSBcYmFye1xib2xkc3ltYm9se1V9fSkuClxlbmR7ZXF1YXRpb259CgoKIyMgQ29uc3RydWN0aW9uIG9mIGFuIGV4YWN0IHNvbHV0aW9uIHRvIHRoZSBvcHRpbWFsIGNvbnRyb2wgcHJvYmxlbSB7I2V4YWN0X3NvbF9mb3Jfb3B0aW1fY29udHJvbF9wcm9ifQoKV2UgZm9sbG93IHRoZSBhcHByb2FjaCBpbiBAR2x1c2EyMDIxZXJyb3Jlc3RpbWF0ZXMgdG8gY29uc3RydWN0IGFuIGV4YWN0IHNvbHV0aW9uIHRvIHRoZSBvcHRpbWFsIGNvbnRyb2wgcHJvYmxlbSBcZXFyZWZ7ZXE6Y29zdGZ1bn0tXGVxcmVme2NvbnRyb2xfY29uc3RyYWludHN9LiBUaGUgc29sdXRpb24gdG8gdGhlIGVsbGlwdGljIHByb2JsZW0KXGJlZ2lue2FsaWdufQogICAgXGxhYmVse2VxOmVsbGlwdGljX3Byb2JsZW19CiAgICBcdGFnezEyfQogICAoXGthcHBhXjIgLSBcRGVsdGFfXEdhbW1hKV57XGFscGhhLzJ9IHUocyk9ZihzKSxccXVhZCBzXGluXEdhbW1hLApcZW5ke2FsaWdufQp3aGVyZSAkdShcY2RvdCkkIHNhdGlzZmllcyB0aGUgS2lyY2hob2ZmIHZlcnRleCBjb25kaXRpb25zIFxlcXJlZntlcTpLY29uZH0sIGlzIGdpdmVuIGJ5ClxiZWdpbnthbGlnbip9CiAgICB1KHMpID0gXHN1bV97aj0xfV57XGluZnR5fVxsYW1iZGFfal57LVxhbHBoYS8yfShmLGVfailfe0xfMihcR2FtbWEpfWVfaihzKSxccXVhZCBzXGluXEdhbW1hLgpcZW5ke2FsaWduKn0KV2UgY2hvb3NlIGNvZWZmaWNpZW50cyAkXHt4X2pcfV97aj0xfV57XGluZnR5fSQgYW5kICRce3lfalx9X3tqPTF9XntcaW5mdHl9JCBhbmQgc2V0ClxiZWdpbnthbGlnbn0KICAgIFxsYWJlbHtlcTpyaHNmdW5jdGlvbnN9CiAgICBcdGFnezEzfQogICBmKHMpID0gXHN1bV97aj0xfV57XGluZnR5fXhfamVfaihzKVxxdWFkIFx0ZXh0eyBhbmQgfSBccXVhZCBnKHMpID0gXHN1bV97aj0xfV57XGluZnR5fXlfamVfaihzKS4KXGVuZHthbGlnbn0KSGVuY2UgdGhlIGNvcnJlc3BvbmRpbmcgc29sdXRpb25zIGFyZSBnaXZlbiBieQpcYmVnaW57YWxpZ259CiAgICBcbGFiZWx7ZXE6c29sc30KICAgIFx0YWd7MTR9CiAgIHUocykgPSBcc3VtX3tqPTF9XntcaW5mdHl9XGxhbWJkYV9qXnstXGFscGhhLzJ9eF9qZV9qKHMpXHF1YWQgXHRleHR7IGFuZCB9IFxxdWFkIHYocykgPSBcc3VtX3tqPTF9XntcaW5mdHl9XGxhbWJkYV9qXnstXGFscGhhLzJ9eV9qZV9qKHMpLgpcZW5ke2FsaWdufQogV2UgYWxzbyBkZWZpbmUgJFxwc2kodCkgPSBcY29zKHQpJCBhbmQgJFxwaGkodCkgPVxzaW4oVC10KSQuIE9ic2VydmUgdGhhdCAkXHBzaSgwKT0xJCBhbmQgJFxwaGkoVCk9MCQuIFdlIGFsc28gY2hvb3NlICRcbXUgPTAuMSQsICRhKHMsdCkgPSAtMC41JCBhbmQgJGIocyx0KT0wLjUkLiBOb3cgd2Ugc2V0ClxiZWdpbnthbGlnbip9CiAgICBmKHMsdCkgJj0gXHBzaScodCl1KHMpK1xwc2kodClmKHMpIC0gXGJhcnt6fShzLHQpID0gLVxzaW4odClcc3VtX3tqPTF9XntcaW5mdHl9XGxhbWJkYV9qXnstXGFscGhhLzJ9eF9qZV9qKHMpKyBcY29zKHQpXHN1bV97aj0xfV57XGluZnR5fXhfamVfaihzKS0gXG1heFx7LTAuNSxcbWluXHswLjUsIFxzaW4oVC10KVxzdW1fe2o9MX1ee1xpbmZ0eX1cbGFtYmRhX2peey1cYWxwaGEvMn15X2plX2oocylcfVx9XFwgCiAgICB1X2Qocyx0KSAmID0gXHBzaSh0KXUocyktXG11XHBoaScodCl2KHMpK1xtdVxwaGkodClnKHMpID0gXGNvcyh0KVxzdW1fe2o9MX1ee1xpbmZ0eX1cbGFtYmRhX2peey1cYWxwaGEvMn14X2plX2oocykrXG11XGNvcyhULXQpXHN1bV97aj0xfV57XGluZnR5fVxsYW1iZGFfal57LVxhbHBoYS8yfXlfamVfaihzKStcbXVcc2luKFQtdClcc3VtX3tqPTF9XntcaW5mdHl9eV9qZV9qKHMpXFwKICAgIHVfMChzKSAgJiA9IHUocyk9IFxzdW1fe2o9MX1ee1xpbmZ0eX1cbGFtYmRhX2peey1cYWxwaGEvMn14X2plX2oocykKXGVuZHthbGlnbip9ClRoZSBleGFjdCBzb2x1dGlvbiB0byB0aGUgb3B0aW1hbCBjb250cm9sIHByb2JsZW0gaXMgZ2l2ZW4gYnkKXGJlZ2lue2FsaWdufQogICAgXGJhcnt1fShzLHQpJj1ccHNpKHQpdShzKSA9IFxjb3ModClcc3VtX3tqPTF9XntcaW5mdHl9XGxhbWJkYV9qXnstXGFscGhhLzJ9eF9qZV9qKHMpXFwKICAgIFxiYXJ7cH0ocyx0KSY9LVxtdVxwaGkodCl2KHMpPSAtIFxtdVxzaW4oVC10KVxzdW1fe2o9MX1ee1xpbmZ0eX1cbGFtYmRhX2peey1cYWxwaGEvMn15X2plX2oocylcXAogICAgXGJhcnt6fShzLHQpJj1cbWF4XHthKHMsdCksXG1pblx7YihzLHQpLCBccGhpKHQpdihzKVx9XH09IFxtYXhcey0wLjUsXG1pblx7MC41LCBcc2luKFQtdClcc3VtX3tqPTF9XntcaW5mdHl9XGxhbWJkYV9qXnstXGFscGhhLzJ9eV9qZV9qKHMpXH1cfSAKXGVuZHthbGlnbn0KCk9ic2VydmUgdGhhdCB0aGUgY2hvaWNlIG9mICRmJCBkZWNvdXBsZXMgdGhlIHN5c3RlbXMgXGVxcmVme251bWVyaWNhbHNjaGVtZTF9LVxlcXJlZnt0aGVudW1lcmljYWxzY2hlbWUyfS4KCiMjIyBDaGVja2luZyB0aGF0ICRcYmFye3V9JCBhbmQgJFxiYXJ7cH0kIHNhdGlzZnkgXGVxcmVme2VxOm1haW5lcW9wdGltYWx9IGFuZCBcZXFyZWZ7ZXE6YWRqb2ludGVxfQoKQXQgJHQ9MCQsIHdlIGdldApcYmVnaW57YWxpZ24qfQogICAgXGJhcnt1fShzLDApID0gXHBzaSgwKXUocykgPSB1KHMpID0gdV8wKHMpLApcZW5ke2FsaWduKn0Kd2hpbGUgZm9yICQocyx0KSBcaW4gXEdhbW1hIFx0aW1lcyAoMCwgVCkkLCB3ZSBvYnRhaW4KXGJlZ2lue2FsaWduKn0KICAgIFxwYXJ0aWFsX3QgXGJhcnt1fShzLHQpICsoXGthcHBhXjIgLSBcRGVsdGFfXEdhbW1hKV57XGFscGhhLzJ9IFxiYXJ7dX0ocyx0KSA9IFxwc2knKHQpdShzKSArIFxwc2kodCkgKFxrYXBwYV4yIC0gXERlbHRhX1xHYW1tYSlee1xhbHBoYS8yfSB1KHMpID0gXHBzaScodCl1KHMpICsgXHBzaSh0KWYocykgPSBmKHMsdCkgKyBcYmFye3p9KHMsdCkuClxlbmR7YWxpZ24qfQpTaW1pbGFybHksIGF0ICR0PVQkLCB3ZSBnZXQKXGJlZ2lue2FsaWduKn0KICAgIFxiYXJ7cH0ocyxUKSA9IC1cbXVccGhpKFQpdihzKSA9IDAsClxlbmR7YWxpZ24qfQp3aGlsZSBmb3IgJChzLHQpIFxpbiBcR2FtbWEgXHRpbWVzICgwLCBUKSQsIHdlIG9idGFpbgpcYmVnaW57YWxpZ24qfQogICAgLVxwYXJ0aWFsX3QgXGJhcntwfShzLHQpICsoXGthcHBhXjIgLSBcRGVsdGFfXEdhbW1hKV57XGFscGhhLzJ9IFxiYXJ7cH0ocyx0KSA9IFxtdVxwaGknKHQpdihzKSAtXG11IFxwaGkodCkgKFxrYXBwYV4yIC0gXERlbHRhX1xHYW1tYSlee1xhbHBoYS8yfSB2KHMpID0gXG11XHBoaScodCl2KHMpIC1cbXUgXHBoaSh0KWcocykKXGVuZHthbGlnbip9CmFuZApcYmVnaW57YWxpZ24qfQogICAgXGJhcnt1fShzLHQpIC0gdV9kKHMsdCkgPSBccHNpKHQpdShzKSAtIFxwc2kodCl1KHMpICsgXG11XHBoaScodCl2KHMpIC0gXG11XHBoaSh0KWcocykgPSBcbXVccGhpJyh0KXYocykgLVxtdSBccGhpKHQpZyhzKS4KXGVuZHthbGlnbip9ClRoaXMgc2hvd3MgdGhhdCAkXGJhcnt1fSQgYW5kICRcYmFye3B9JCBzYXRpc2Z5IFxlcXJlZntlcTptYWluZXFvcHRpbWFsfSBhbmQgXGVxcmVme2VxOmFkam9pbnRlcX0uCgojIyBOdW1lcmljYWwgaW1wbGVtZW50YXRpb24geyNudW1faW1wbGVtZW50YXRpb25fb3B0aW1hbF9jb250cm9sX3Byb2JsZW19CgpJbiB0aGlzIHNlY3Rpb24sIHdlIGltcGxlbWVudCB0aGUgbnVtZXJpY2FsIHNjaGVtZSBmb3IgdGhlIG9wdGltYWwgY29udHJvbCBwcm9ibGVtIFxlcXJlZntlcTpjb3N0ZnVufS1cZXFyZWZ7Y29udHJvbF9jb25zdHJhaW50c30gb24gdGhlIHRhZHBvbGUgZ3JhcGguIFdlIHN0YXJ0IGJ5IGRlZmluaW5nIHRoZSBwYXJhbWV0ZXJzIGFuZCBjb25zdHJ1Y3RpbmcgdGhlIGV4YWN0IHNvbHV0aW9uLiBXZSB0aGVtIHVzZSBcZXFyZWZ7bnVtZXJpY2Fsc2NoZW1lMX0tXGVxcmVme3RoZW51bWVyaWNhbHNjaGVtZTJ9IHRvIHNvbHZlIHRoZSBmb3J3YXJkIGFuZCBhZGpvaW50IHByb2JsZW1zLiBGaW5hbGx5LCB3ZSBjb21wdXRlIHRoZSBjb250cm9sIHZhcmlhYmxlICRcYmFye3p9JCB1c2luZyBcZXFyZWZ7enp9LgoKCmBgYHtyfQojIFBhcmFtZXRlcnMKVF9maW5hbCA8LSAyCnRpbWVfc3RlcCA8LSAwLjAyCmggPC0gMC4wMgprYXBwYSA8LSAxNQphbHBoYSA8LSAxLjgKbSA9IDQKYmV0YSA8LSBhbHBoYS8yCm11IDwtIDAuMQphIDwtIC0gMC41CmIgPC0gMC41Ck5fZmluaXRlID0gNCAjIGNob29zZSBldmVuCmFkanVzdGVkX05fZmluaXRlIDwtIE5fZmluaXRlICsgTl9maW5pdGUvMiArIDEKIyBDb2VmZmljaWVudHMgZm9yIGYgYW5kIGcKY29lZmZfZWxsaXB0aWNfZyA8LSByZXAoMCwgYWRqdXN0ZWRfTl9maW5pdGUpCmNvZWZmX2VsbGlwdGljX2dbN10gPC0gMTAwCmNvZWZmX2VsbGlwdGljX2YgPC0gcmVwKDAsIGFkanVzdGVkX05fZmluaXRlKQpjb2VmZl9lbGxpcHRpY19mWzddIDwtIC0xMAojIHRpbWUgYW5kIHNwYXRpYWwgbWVzaCAKdGltZV9zZXEgPC0gc2VxKDAsIFRfZmluYWwsIGxlbmd0aC5vdXQgPSAoKFRfZmluYWwgLSAwKSAvIHRpbWVfc3RlcCArIDEpKQpncmFwaCA8LSBnZXRzLmdyYXBoLnRhZHBvbGUoaCA9IGgpCmdyYXBoJGNvbXB1dGVfZmVtKCkKRyA8LSBncmFwaCRtZXNoJEcKQyA8LSBncmFwaCRtZXNoJEMKTCA8LSBrYXBwYV4yKkMgKyBHCkkgPC0gTWF0cml4OjpEaWFnb25hbChucm93KEMpKQojIENvbXB1dGUgdGhlIGVpZ2VudmFsdWVzIGFuZCBlaWdlbmZ1bmN0aW9ucwplaWdlbl9wYXJhbXMgPC0gZ2V0cy5laWdlbi5wYXJhbXMoTl9maW5pdGUgPSBOX2Zpbml0ZSwgCiAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICBrYXBwYSA9IGthcHBhLCAKICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgIGFscGhhID0gYWxwaGEsIAogICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgZ3JhcGggPSBncmFwaCkKRUlHRU5WQUxfTUlOVVNfQUxQSEEgPC0gZWlnZW5fcGFyYW1zJEVJR0VOVkFMX01JTlVTX0FMUEhBCkVJR0VORlVOIDwtIGVpZ2VuX3BhcmFtcyRFSUdFTkZVTgojIENvbnN0cnVjdCB0aGUgcmlnaHQgaGFuZCBzaWRlIGZ1bmN0aW9ucyBmIGFuZCBnIGZvciB0aGUgZWxsaXB0aWMgcHJvYmxlbQplbGxpcHRpY19mIDwtIGFzLnZlY3RvcihFSUdFTkZVTiAlKiUgY29lZmZfZWxsaXB0aWNfZikKZWxsaXB0aWNfZyA8LSBhcy52ZWN0b3IoRUlHRU5GVU4gJSolIGNvZWZmX2VsbGlwdGljX2cpCiMgQ29uc3RydWN0IHRoZSBjb3JyZXNwb25kaW5nIGVsbGlwdGljIHNvbHV0aW9uIHUgYW5kIHYKZWxsaXB0aWNfdSA8LSBhcy52ZWN0b3IoRUlHRU5GVU4gJSolIChjb2VmZl9lbGxpcHRpY19mICogRUlHRU5WQUxfTUlOVVNfQUxQSEEpKQplbGxpcHRpY192IDwtIGFzLnZlY3RvcihFSUdFTkZVTiAlKiUgKGNvZWZmX2VsbGlwdGljX2cgKiBFSUdFTlZBTF9NSU5VU19BTFBIQSkpCgpBIDwtIG1hdHJpeChhLCBucm93ID0gbnJvdyhDKSwgbmNvbCA9IGxlbmd0aCh0aW1lX3NlcSkpCkIgPC0gbWF0cml4KGIsIG5yb3cgPSBucm93KEMpLCBuY29sID0gbGVuZ3RoKHRpbWVfc2VxKSkKIyBDb25zdHJ1Y3QgdGhlIG9wdGltYWwgdmFyaWFibGVzCnBzaSA8LSBjb3ModGltZV9zZXEpCnBzaV9wcmltZSA8LSAtIHNpbih0aW1lX3NlcSkKcGhpIDwtIHNpbihUX2ZpbmFsIC0gdGltZV9zZXEpCnBoaV9wcmltZSA8LSAtIGNvcyhUX2ZpbmFsIC0gdGltZV9zZXEpCiMgcHNpX3JhdGUgPC0gMQojIHBoaV9yYXRlIDwtIDEKIyBwaGlfY29lZmZpY2llbnQgPC0gMQojIHBzaSA8LSBleHAocHNpX3JhdGUqdGltZV9zZXEpCiMgcHNpX3ByaW1lIDwtIHBzaV9yYXRlICogZXhwKHBzaV9yYXRlKnRpbWVfc2VxKQojIHBoaSA8LSBwaGlfY29lZmZpY2llbnQgKiBleHAocGhpX3JhdGUqdGltZV9zZXEpIC0gcGhpX2NvZWZmaWNpZW50ICogZXhwKHBoaV9yYXRlKlRfZmluYWwpCiMgcGhpX3ByaW1lIDwtIHBoaV9jb2VmZmljaWVudCAqIHBoaV9yYXRlICogZXhwKHBoaV9yYXRlKnRpbWVfc2VxKQp1X2JhciA8LSBlbGxpcHRpY191ICUqJSB0KHBzaSkKcF9iYXIgPC0gLSBtdSAqIChlbGxpcHRpY192ICUqJSB0KHBoaSkpCnpfYmFyIDwtIHBtYXgoQSwgcG1pbihCLCAtIHBfYmFyIC8gbXUpKQojIENvbnN0cnVjdCB0aGUgZGF0YSBmb3IgdGhlIHByb2JsZW0KZiA8LSBlbGxpcHRpY191ICUqJSB0KHBzaV9wcmltZSkgKyBlbGxpcHRpY19mICUqJSB0KHBzaSkgLSB6X2JhcgpmX3BsdXNfel9iYXIgPC0gZWxsaXB0aWNfdSAlKiUgdChwc2lfcHJpbWUpICsgZWxsaXB0aWNfZiAlKiUgdChwc2kpCnVfZCA8LSBlbGxpcHRpY191ICUqJSB0KHBzaSkgLQogIG11ICogKGVsbGlwdGljX3YgJSolIHQocGhpX3ByaW1lKSkgKwogIG11ICogKGVsbGlwdGljX2cgJSolIHQocGhpKSkKdV8wIDwtIGVsbGlwdGljX3UKIyBDb25zdHJ1Y3QgdGhlIGZyYWN0aW9uYWwgb3BlcmF0b3IsIHdoaWNoIGlzIHNoYXJlZCBmb3IgdGhlIGZvcndhcmQgYW5kIGFkam9pbnQgcHJvYmxlbXMKbXlfb3BfZnJhYyA8LSBteS5mcmFjdGlvbmFsLm9wZXJhdG9ycy5mcmFjKEwsIGJldGEsIEMsIHNjYWxlLmZhY3RvciA9IGthcHBhXjIsIG0gPSBtLCB0aW1lX3N0ZXApCiMgQ29uc3RydWN0IHRoZSBpbnRlZ3JhdGlvbiBtZXNoCm92ZXJraWxsX2dyYXBoIDwtIGdldHMuZ3JhcGgudGFkcG9sZShoID0gMC4wMDAxKQpvdmVya2lsbF9ncmFwaCRjb21wdXRlX2ZlbSgpCm92ZXJraWxsX0VJR0VORlVOIDwtIGdldHMuZWlnZW4ucGFyYW1zKE5fZmluaXRlID0gTl9maW5pdGUsIAogICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICBrYXBwYSA9IGthcHBhLCAKICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgYWxwaGEgPSBhbHBoYSwgCiAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgIGdyYXBoID0gb3ZlcmtpbGxfZ3JhcGgpJEVJR0VORlVOCiMgQ29uc3RydWN0IHRoZSByaWdodCBoYW5kIHNpZGUgZnVuY3Rpb25zIGYgYW5kIGcgZm9yIHRoZSBlbGxpcHRpYyBwcm9ibGVtIG9uIHRoZSBpbnRlZ3JhdGlvbiBtZXNoCm92ZXJraWxsX2VsbGlwdGljX2YgPC0gYXMudmVjdG9yKG92ZXJraWxsX0VJR0VORlVOICUqJSBjb2VmZl9lbGxpcHRpY19mKQpvdmVya2lsbF9lbGxpcHRpY19nIDwtIGFzLnZlY3RvcihvdmVya2lsbF9FSUdFTkZVTiAlKiUgY29lZmZfZWxsaXB0aWNfZykKIyBDb25zdHJ1Y3QgdGhlIGNvcnJlc3BvbmRpbmcgZWxsaXB0aWMgc29sdXRpb24gdSBhbmQgdiBvbiB0aGUgaW50ZWdyYXRpb24gbWVzaApvdmVya2lsbF9lbGxpcHRpY191IDwtIGFzLnZlY3RvcihvdmVya2lsbF9FSUdFTkZVTiAlKiUgKGNvZWZmX2VsbGlwdGljX2YgKiBFSUdFTlZBTF9NSU5VU19BTFBIQSkpCm92ZXJraWxsX2VsbGlwdGljX3YgPC0gYXMudmVjdG9yKG92ZXJraWxsX0VJR0VORlVOICUqJSAoY29lZmZfZWxsaXB0aWNfZyAqIEVJR0VOVkFMX01JTlVTX0FMUEhBKSkKIyBDb25zdHJ1Y3QgdGhlIHByb2plY3Rpb24gbWF0cml4ClBzaSA8LSBncmFwaCRmZW1fYmFzaXMob3ZlcmtpbGxfZ3JhcGgkZ2V0X21lc2hfbG9jYXRpb25zKCkpClIgPC0gdChQc2kpICUqJSBvdmVya2lsbF9ncmFwaCRtZXNoJEMKIyBDb25zdHJ1Y3QgdGhlIHJpZ2h0IGhhbmQgc2lkZSB0ZXJtIGZvciB0aGUgZm9yd2FyZCBwcm9ibGVtCkZfcGx1c19aX2JhciA8LSBvdmVya2lsbF9lbGxpcHRpY191ICUqJSB0KHBzaV9wcmltZSkgKyBvdmVya2lsbF9lbGxpcHRpY19mICUqJSB0KHBzaSkgCkZfcGx1c19aX2Jhcl9pbm5lcnByb2R1Y3QgPC0gUiAlKiUgRl9wbHVzX1pfYmFyCiMgU29sdmUgdGhlIGZvcndhcmQgcHJvYmxlbQpVX2JhciA8LSBzb2x2ZV9mcmFjdGlvbmFsX2V2b2x1dGlvbihteV9vcF9mcmFjLCAKICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgIHRpbWVfc3RlcCwgCiAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICB0aW1lX3NlcSwgCiAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICB2YWxfYXRfMCA9IHVfMCwgCiAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICBSSFNUID0gRl9wbHVzX1pfYmFyX2lubmVycHJvZHVjdCkKIyBDb25zdHJ1Y3QgdGhlIHJpZ2h0IGhhbmQgc2lkZSB0ZXJtIGZvciB0aGUgYWRqb2ludCBwcm9ibGVtClVfZCA8LSBvdmVya2lsbF9lbGxpcHRpY191ICUqJSB0KHBzaSkgLQogIG11ICogKG92ZXJraWxsX2VsbGlwdGljX3YgJSolIHQocGhpX3ByaW1lKSkgKyAKICBtdSAqIChvdmVya2lsbF9lbGxpcHRpY19nICUqJSB0KHBoaSkpClZfZCA8LSByZXZlcnNlY29sdW1ucyhSICUqJSBVX2QpClZfYmFyIDwtIHJldmVyc2Vjb2x1bW5zKFIgJSolIFBzaSAlKiUgVV9iYXIpCiMgU29sdmUgdGhlIGFkam9pbnQgcHJvYmxlbQpRX2JhciA8LSBzb2x2ZV9mcmFjdGlvbmFsX2V2b2x1dGlvbihteV9vcF9mcmFjLCAKICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgIHRpbWVfc3RlcCwgCiAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICB0aW1lX3NlcSwgCiAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICB2YWxfYXRfMCA9IHVfMCowLCAKICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgICAgIFJIU1QgPSBWX2JhciAtIFZfZCkKUF9iYXIgPC0gcmV2ZXJzZWNvbHVtbnMoUV9iYXIpCiMgQ29tcHV0ZSB0aGUgY29udHJvbCB2YXJpYWJsZQpaX2JhciA8LSBwbWF4KEEsIHBtaW4oQiwgLSBQX2JhciAvIG11KSkKYGBgCgoKYGBge3J9CnN0YXJ0X2lkeCA8LSB3aGljaC5taW4oYWJzKHRpbWVfc2VxIC0gMCkpCmVuZF9pZHggPC0gd2hpY2gubWluKGFicyh0aW1lX3NlcSAtIFRfZmluYWwpKQppZHggPC0gc2VxKHN0YXJ0X2lkeCwgZW5kX2lkeCwgYnkgPSAxKQpgYGAKCmBgYHtyLCBmaWcuaGVpZ2h0ID0gNiwgb3V0LndpZHRoID0gIjEwMCUiLCBmaWcuY2FwID0gY2FwdGlvbmVyKCJSaWdodC1oYW5kIHNpZGUgZnVuY3Rpb25zICRmJCBhbmQgJGckIGZvciB0aGUgZWxsaXB0aWMgcHJvYmxlbS4iKX0KZ3JhcGgucGxvdHRlci4zZC5zdGF0aWMoZ3JhcGgsIGVsbGlwdGljX2YsIGVsbGlwdGljX2cpCmBgYAoKCmBgYHtyLCBmaWcuaGVpZ2h0ID0gNiwgb3V0LndpZHRoID0gIjEwMCUiLCBmaWcuY2FwID0gY2FwdGlvbmVyKCJTb2x1dGlvbnMgJHUkIGFuZCAkdiQgdG8gdGhlIGVsbGlwdGljIHByb2JsZW0uIil9CmdyYXBoLnBsb3R0ZXIuM2Quc3RhdGljKGdyYXBoLCBlbGxpcHRpY191LCBlbGxpcHRpY192KQpgYGAKCmBgYHtyLCBmaWcuaGVpZ2h0ID0gNiwgb3V0LndpZHRoID0gIjEwMCUiLCBmaWcuY2FwID0gY2FwdGlvbmVyKCJUaW1lIGV2b2x1dGlvbiBvZiB0aGUgcmlnaHQtaGFuZCBzaWRlIGZ1bmN0aW9uICRmK1xcYmFye3p9JCBvZiB0aGUgc3RhdGUgZXF1YXRpb24uIil9CmdyYXBoLnBsb3R0ZXIuM2Quc2luZ2xlKGdyYXBoLCBmX3BsdXNfel9iYXJbLCBpZHhdLCB0aW1lX3NlcVtpZHhdKQpgYGAKCmBgYHtyLCBmaWcuaGVpZ2h0ID0gNiwgb3V0LndpZHRoID0gIjEwMCUiLCBmaWcuY2FwID0gY2FwdGlvbmVyKCJUaW1lIGV2b2x1dGlvbiBvZiB0aGUgcmlnaHQtaGFuZCBzaWRlIGZ1bmN0aW9uICRcXGJhcnt2fS12X2QkIG9mIHRoZSBhZGpvaW50IGVxdWF0aW9uLiIpfQp2X2Jhcl9taW51c192X2QgPC0gcmV2ZXJzZWNvbHVtbnModV9iYXIgLSB1X2QpCmdyYXBoLnBsb3R0ZXIuM2Quc2luZ2xlKGdyYXBoLCB2X2Jhcl9taW51c192X2RbLCBpZHhdLCB0aW1lX3NlcVtpZHhdKQpgYGAKCiMjIFZpc3VhbGl6YXRpb24gb2YgdGhlIHJlc3VsdHMgeyN2aXN1YWxpemF0aW9uX3Jlc3VsdHN9CgpgYGB7ciwgZmlnLmhlaWdodCA9IDYsIG91dC53aWR0aCA9ICIxMDAlIiwgZmlnLmNhcCA9IGNhcHRpb25lcigiRXhhY3QgKGluIHJlZCkgYW5kIG51bWVyaWNhbCAoaW4gYmx1ZSkgb3B0aW1hbCBzdGF0ZSB2YXJpYWJsZS4iKX0KZ3JhcGgucGxvdHRlci4zZChncmFwaCwgdGltZV9zZXFbaWR4XSwgdGltZV9zZXFbaWR4XSwgdV9iYXJbLCBpZHhdLCBVX2JhclssIGlkeF0pCmBgYAoKCgpgYGB7ciwgZmlnLmhlaWdodCA9IDYsIG91dC53aWR0aCA9ICIxMDAlIiwgZmlnLmNhcCA9IGNhcHRpb25lcigiRXhhY3QgKGluIHJlZCkgYW5kIG51bWVyaWNhbCAoaW4gYmx1ZSkgb3B0aW1hbCBhZGpvaW50IHZhcmlhYmxlLiIpfQpncmFwaC5wbG90dGVyLjNkKGdyYXBoLCB0aW1lX3NlcVtpZHhdLCB0aW1lX3NlcVtpZHhdLCBwX2JhclssIGlkeF0sIFBfYmFyWywgaWR4XSkKYGBgCgoKYGBge3IsIGZpZy5oZWlnaHQgPSA2LCBvdXQud2lkdGggPSAiMTAwJSIsIGZpZy5jYXAgPSBjYXB0aW9uZXIoIkV4YWN0IChpbiByZWQpIGFuZCBudW1lcmljYWwgKGluIGJsdWUpIG9wdGltYWwgY29udHJvbCB2YXJpYWJsZS4iKX0KZ3JhcGgucGxvdHRlci4zZChncmFwaCwgdGltZV9zZXFbaWR4XSwgdGltZV9zZXFbaWR4XSwgel9iYXJbLCBpZHhdLCBaX2JhclssIGlkeF0pCmBgYAoKCmBgYHtyLCBmaWcuaGVpZ2h0ID0gNiwgb3V0LndpZHRoID0gIjEwMCUiLCBmaWcuY2FwID0gY2FwdGlvbmVyKCJEZXNpcmVkIHN0YXRlIChpbiByZWQpIGFuZCBvcHRpbWFsIHN0YXRlIChpbiBibHVlKS4iKX0KZ3JhcGgucGxvdHRlci4zZChncmFwaCwgdGltZV9zZXFbaWR4XSwgdGltZV9zZXFbaWR4XSwgdV9kWywgaWR4XSwgdV9iYXJbLCBpZHhdKQpgYGAKCgojIFJlZmVyZW5jZXMKCmBgYHtyfQpjaXRlX3BhY2thZ2VzKG91dHB1dCA9ICJwYXJhZ3JhcGgiLCBvdXQuZGlyID0gIi4iKQpgYGA=