TransportMaps.Distributions.Examples.BurgersPDE.BurgersProblems

Module Contents

Classes

BurgersSolver

Defines the solver for the Burger's problem with respect to prescribed initial and boundary conditions and viscosity

AdjointFinalL2MismatchBurgersSolver

Defines the adjoint for the Burger's problem with respect to initial conditions and viscosity

FinalL2MismatchBurgersProblem

Defines the adjoint for the Burger's problem with respect to initial conditions and viscosity

class TransportMaps.Distributions.Examples.BurgersPDE.BurgersProblems.BurgersSolver(**kwargs)[source]

Bases: TransportMaps.DOLFIN.Solver

Defines the solver for the Burger’s problem with respect to prescribed initial and boundary conditions and viscosity

It discretizes the time dependent problem using the backward Euler method.

\[\begin{split}\begin{cases} \partial_t u({\bf x},t) = \nabla \cdot (\mu \nabla u({\bf x},t)) - u({\bf x},t) \sum_{i=1}^d \partial_{x_i}u({\bf x},t) & {\bf x} \in \Omega \\ u({\bf x},0) = u_0({\bf x}) & t = 0\\ u({\bf x},t) = g({\bf x}) & {\bf x} \in \partial\Omega \end{cases}\end{split}\]
set_up(**kwargs)[source]
Kwargs:

VEFS (FunctionSpace): function space where the problem is defined on dt (float): time step bcs (list): list of boundary conditions nsteps (int): number of integration steps

solve(nu, u0)[source]
Parameters:
  • nu (float) – viscosity

  • u0 (ndarray) – initial conditions

Returns:

(ndarray) – solution \(u({\bf x}, T)\)

class TransportMaps.Distributions.Examples.BurgersPDE.BurgersProblems.AdjointFinalL2MismatchBurgersSolver(**kwargs)[source]

Bases: BurgersSolver

Defines the adjoint for the Burger’s problem with respect to initial conditions and viscosity

It discretizes the time dependent problem using the backward Euler method.

\[\begin{split}\begin{cases} \partial_t u({\bf x},t) = \nabla \cdot (\mu \nabla u({\bf x},t)) - u({\bf x},t) \sum_{i=1}^d \partial_{x_i}u({\bf x},t) & {\bf x} \in \Omega \\ u({\bf x},0) = u_0({\bf x}) & t = 0\\ u({\bf x},t) = g({\bf x}) & {\bf x} \in \partial\Omega \end{cases}\end{split}\]

and compute the functional \(J(u) = \int (u(T) - u_d)^2 dx\), along with the gradients \(\partial_{u_0} J\) and \(\partial_\nu J\).

property ud[source]
set_up(ud=None, **kwargs)[source]

If the optional keyword arguments are provided, then builds the reduced functional to the be used in function evaluations. To do so, an forward solve must be run, so that doflin_adjoint can gather all the list of operations.

Optional Kwargs:

ud (ndarray): data forthe final tie \(u_d\)

new_function()[source]

Overrides the function generator to created dolfin_adjoint.Function instead of ordinary dolfin.Function.

__getstate__()[source]

Avoids storing the dolfin_adjoint.ReducedFunctional Jhat

__setstate__(state)[source]

Reloads the dolfin_adjoint.ReducedFunctional Jhat

_set_up_reduced_functional()[source]
evaluate(nu, u0)[source]
Parameters:
  • nu (float) – viscosity \(\nu\)

  • u0 (ndarray) – initial conditions \(u_0\)

Returns

(float): \(J\)

grad_x(nu, u0)[source]
Parameters:
  • nu (float) – viscosity \(\nu\)

  • u0 (ndarray) – initial conditions \(u_0\)

Returns:

(tuple) – \(( \partial_{u_0} J, \partial_{\nu} J )\)

Note

since the adjoint is used, the forward model evaluation must also be called. Therefore this function behaves exactly as :fun:`tuple_grad_x`, but returns only the gradient.

tuple_grad_x(nu, u0)[source]
Parameters:
  • nu (float) – viscosity \(\nu\)

  • u0 (ndarray) – initial conditions \(u_0\)

Returns:

(tuple) – \(( J, \partial_{\nu} J, \partial_{u_0} J )\)

class TransportMaps.Distributions.Examples.BurgersPDE.BurgersProblems.FinalL2MismatchBurgersProblem(**kwargs)[source]

Bases: AdjointFinalL2MismatchBurgersSolver

Defines the adjoint for the Burger’s problem with respect to initial conditions and viscosity

It discretizes the time dependent problem using the backward Euler method.

\[\begin{split}\begin{cases} \partial_t u({\bf x},t) = \nabla \cdot (\mu \nabla u({\bf x},t)) - u({\bf x},t) \sum_{i=1}^d \partial_{x_i}u({\bf x},t) & {\bf x} \in \Omega \\ u({\bf x},0) = u_0({\bf x}) & t = 0\\ u({\bf x},t) = g({\bf x}) & {\bf x} \in \partial\Omega \end{cases}\end{split}\]

and compute the functional \(J(u) = \int (u(T) - u_d)^2 dx\), along with the gradients \(\partial_{u_0} J\) and \(\partial_\nu J\).

_init_args = ['T', 'xl', 'xr', 'ul', 'ur', 'nels', 'order', 'nsteps'][source]
__getstate__()[source]

Avoids storing the dolfin_adjoint.ReducedFunctional Jhat

set_up()[source]

If the optional keyword arguments are provided, then builds the reduced functional to the be used in function evaluations. To do so, an forward solve must be run, so that doflin_adjoint can gather all the list of operations.

Optional Kwargs:

ud (ndarray): data forthe final tie \(u_d\)