soupy.modeling

soupy.modeling.PDEControlProblem

class soupy.modeling.PDEControlProblem.PDEVariationalControlProblem(*args: Any, **kwargs: Any)[source]

Bases: PDEVariationalProblem

Constructor

Parameters:
  • Vh (list of dolfin.FunctionSpace) – List of function spaces the state, parameter, adjoint, and control

  • varf_handler – Variational form handler with __call__ method

  • bc – List of Dirichlet boundary conditions for the state

  • bc0 – List of zeroed Dirichlet boundary conditions

  • is_fwd_linear (bool) – Flag indicating whether the forward problem is linear

  • lu_method (str) – Method for solving linear systems (default, mumps, etc.)

apply_ij(i, j, dir, out)[source]

Given \(u, m, p, z\); compute \(\delta_{ij} F(u, m, p, z; \hat{i}, \tilde{j})\) in the direction \(\tilde{j} =\) dir, \(\forall \hat{i}\).

apply_ijk(i, j, k, x, jdir, kdir, out)[source]
evalGradientControl(x, out)[source]

Given \(u, m, p, z\); evaluate \(\delta_z F(u, m, p, z; \hat{z}),\, \forall \hat{z}.\)

evalGradientParameter(x, out)[source]

Given \(u, m, p, z\); evaluate \(\delta_m F(u, m, p, z; \hat{m}),\, \forall \hat{m}.\)

generate_control()[source]
generate_parameter()[source]

Return a vector in the shape of the parameter.

generate_state()[source]

Return a vector in the shape of the state.

init_control(z)[source]

Initialize the parameter.

init_parameter(m)[source]

Initialize the parameter.

setLinearizationPoint(x, gauss_newton_approx)[source]

Set the values of the state and parameter for the incremental forward and adjoint solvers.

set_nonlinear_solver_parameters(parameters)[source]

Set the solver parameters used for dolfin.NonlinearVariationalSolver

Parameters:

parameters (dict) – Solver parameters for dolfin.NonlinearVariationalSolver

solveAdj(adj, x, adj_rhs)[source]

Solve the linear adjoint problem: Given \(m, z, u\); find \(p\) such that

\[\delta_u F(u, m, p, z;\hat{u}) = 0, \quad \forall \hat{u}.\]
solveFwd(state, x)[source]

Solve the possibly nonlinear forward problem: Given \(m, z\), find \(u\) such that

\[\delta_p F(u, m, p, z;\hat{p}) = 0,\quad \forall \hat{p}.\]
solveIncremental(out, rhs, is_adj)[source]

If is_adj == False:

Solve the forward incremental system: Given \(u, m, z\), find \(\tilde{u}\) such that

\[\delta_{pu} F(u, m, p, z ; \hat{p}, \tilde{u}) = \mbox{rhs},\quad \forall \hat{p}.\]

If is_adj == True:

Solve the adjoint incremental system: Given \(u, m, z\), find \(\tilde{p}\) such that

\[\delta_{up} F(u, m, p, z; \hat{u}, \tilde{p}) = \mbox{rhs},\quad \forall \hat{u}.\]

soupy.modeling.controlModel

class soupy.modeling.controlModel.ControlModel(problem, qoi)[source]

Bases: object

This class contains the structure needed to evaluate the control objective As inputs it takes a PDEVariationalControlProblem, and a Qoi object.

In the following we will denote with
  • u the state variable

  • m the (model) parameter variable

  • p the adjoint variable

  • z the control variable

Constructor

Parameters:
  • problem (PDEVariationalControlProblem) – The PDE problem

  • qoi (ControlQoI) – The quantity of interest

applyC(dm, out)[source]

Apply the \(C_{m}\) block of the Hessian to a (incremental) parameter variable, i.e. out = \(C_{m} dm\)

Parameters:
  • dm (dolfin.Vector) – The (incremental) parameter variable

  • out (dolfin.Vector) – The action of the \(C_z\) block on dm

Note

This routine assumes that out has the correct shape.

applyCt(dp, out)[source]

Apply the transpose of the \(C_{m}\) block of the Hessian to a (incremental) adjoint variable. out = \(C_{m}^T dp\)

Parameters:
  • dp (dolfin.Vector) – The (incremental) adjoint variable

  • out (dolfin.Vector) – The action of the \(C_{m}^T\) block on dp

..note:: This routine assumes that out has the correct shape.

applyCz(dz, out)[source]

Apply the \(C_z\) block of the Hessian to a (incremental) control variable, i.e. out = \(C_z dz\)

Parameters:
  • dz (dolfin.Vector) – The (incremental) control variable

  • out (dolfin.Vector) – The action of the \(C_z\) block on dz

Note

This routine assumes that out has the correct shape.

applyCzt(dp, out)[source]

Apply the transpose of the \(C_{z}\) block of the Hessian to a (incremental) adjoint variable. out = \(C_{z}^T dp\)

Parameters:
  • dp (dolfin.Vector) – The (incremental) adjoint variable

  • out (dolfin.Vector) – The action of the \(C_{z}^T\) block on dp

..note:: This routine assumes that out has the correct shape.

applyWmm(dm, out)[source]

Apply the \(W_{mm}\) block of the Hessian to a (incremental) parameter variable. out = \(W_{mm} dm\)

Parameters:
  • dm (dolfin.Vector) – The (incremental) parameter variable

  • out (dolfin.Vector) – The action of the \(W_{mm}\) block on du

Note

This routine assumes that out has the correct shape.

applyWmu(du, out)[source]

Apply the \(W_{mu}\) block of the Hessian to a (incremental) state variable. out = \(W_{mu} du\)

Parameters:
  • du (dolfin.Vector) – The (incremental) state variable

  • out (dolfin.Vector) – The action of the \(W_{mu}\) block on du

Note

This routine assumes that out has the correct shape.

applyWum(dm, out)[source]

Apply the \(W_{um}\) block of the Hessian to a (incremental) parameter variable. out = \(W_{um} dm\)

Parameters:
  • dm (dolfin.Vector) – The (incremental) parameter variable

  • out (dolfin.Vector) – The action of the \(W_{um}\) block on du

Note

This routine assumes that out has the correct shape.

applyWuu(du, out)[source]

Apply the \(W_{uu}\) block of the Hessian to a (incremental) state variable. out = \(W_{uu} du\)

Parameters:
  • du (dolfin.Vector) – The (incremental) state variable

  • out (dolfin.Vector) – The action of the \(W_{uu}\) block on du

Note

This routine assumes that out has the correct shape.

applyWuz(dz, out)[source]

Apply the \(W_{uz}\) block of the Hessian to a (incremental) control variable. out = \(W_{uz} dz\)

Parameters:
  • dz (dolfin.Vector) – The (incremental) control variable

  • out (dolfin.Vector) – The action of the \(W_{uz}\) block on du

Note

This routine assumes that out has the correct shape.

applyWzu(du, out)[source]

Apply the \(W_{zu}\) block of the Hessian to a (incremental) state variable. out = \(W_{zu} du\)

Parameters:
  • du (dolfin.Vector) – The (incremental) state variable

  • out (dolfin.Vector) – The action of the \(W_{zu}\) block on du

Note

This routine assumes that out has the correct shape.

applyWzz(dz, out)[source]

Apply the \(W_{zz}\) block of the Hessian to a (incremental) control variable. out = \(W_{zz} dz\)

Parameters:
  • dz (dolfin.Vector) – The (incremental) control variable

  • out (dolfin.Vector) – The action of the \(W_{zz}\) block on du

Note

This routine assumes that out has the correct shape.

apply_ij(i, j, d, out)[source]

Apply the \((i,j)\) block of the Hessian to a vector d

Parameters:
  • i (int) – The output variable index soupy.STATE, soupy.PARAMETER, soupy.ADJOINT, or soupy.CONTROL

  • j (int) – The input variable index soupy.STATE, soupy.PARAMETER, soupy.ADJOINT, or soupy.CONTROL

  • d (dolfin.Vector) – The vector to which the Hessian is applied

  • out (dolfin.Vector) – The action of the Hessian on d

cost(x)[source]

Evaluate the QoI at the a given point \(q(u,m,z)\).

Parameters:

x (list of dolfin.Vector objects) – The point x = [u,m,p,z] at which to evaluate the QoI

Returns:

The value of the QoI

Note

p is not needed to compute the cost functional

evalGradientControl(x, mg)[source]

Evaluate the \(z\) gradient action form at the point \((u,m,p,z)\)

Parameters:
  • x (list of dolfin.Vector objects) – The point x = [u,m,p,z]. Provides the state, parameter, adjoint, and control variables u, m, p, and z for the assembling the gradient action Vector u is also used to assemble the adjoint right hand side.

  • mg (dolfin.Vector) – Dual of the gradient with respect to the control i.e. \((g_z, z_{\mathrm{test}})_{Z}\)

Returns:

the norm of the gradient in the correct inner product \((g_z,g_z)_{Z}^{1/2}\)

evalGradientParameter(x, mg)[source]

Evaluate the \(m\) gradient action form at the point \((u,m,p,z)\)

Parameters:
  • x (list of dolfin.Vector objects) – The point x = [u,m,p,z]. Provides the state, parameter, adjoint, and control variables u, m, p, and z for the assembling the gradient action. Vector u is also used to assemble the adjoint right hand side.

  • mg (dolfin.Vector) – Dual of the gradient with respect to the parameter i.e. \((g_m, m_{\mathrm{test}})_{M}\)

Returns:

the norm of the gradient in the correct inner product \((g_m,g_m)_{M}^{1/2}\)

generate_vector(component='ALL')[source]
Parameters:

component – The component of the vector to generate (soupy.STATE, soupy.PARAMETER, soupy.ADJOINT, soupy.CONTROL, or "ALL")

Returns:

By default, component == "ALL" will return the list [u,m,p,z] where each element is a dolfin.Vector object of the appropriate size If component = soupy.STATE return only u. If component = soupy.PARAMETER return only m. If component = soupy.ADJOINT return only p. If component = soupy.CONTROL return only z.

init_control(z)[source]

Reshape z so that it is compatible with the control variable

init_parameter(m)[source]

Reshape m so that it is compatible with the parameter variable

setLinearizationPoint(x, gauss_newton_approx=False)[source]

Specify the point x = [u,m,p,z] at which the Hessian operator (or the Gauss-Newton approximation) needs to be evaluated.

Parameters:
  • x – The point x = [u,m,p,z] for which the Hessian needs to be evaluated

  • gauss_newton_approx (bool) – whether to use the Gauss-Newton approximation (default: use Newton)

Note

This routine should either: 1. simply store a copy of x and evaluate action of blocks of the Hessian on the fly, or 2. partially precompute the block of the hessian (if feasible)

solveAdj(out, x)[source]

Solve the linear adjoint problem.

Parameters:
  • out (dolfin.Vector) – Solution of the forward problem (state)

  • x – The point x = [u,m,p,z]. Provides the state, parameter and control variables u, m, and z for the solution assembling the adjoint operator. Vector u is also used to assemble the adjoint right hand side.

Type:

list of dolfin.Vector objects

Note

p is not accessed

solveAdjIncremental(sol, rhs)[source]

Solve the incremental adjoint problem for a given right-hand side

Parameters:
  • sol (dolfin.Vector) – Solution of the incremental adjoint problem (adjoint)

  • rhs (dolfin.Vector) – Right hand side of the linear system

solveFwd(out, x)[source]

Solve the (possibly non-linear) forward problem.

Parameters:
  • out (dolfin.Vector) – Solution of the forward problem (state)

  • x – The point x = [u,m,p,z]. Provides the parameter and control variables m and z for the solution of the forward problem and the initial guess u if the forward problem is non-linear

Type:

list of dolfin.Vector objects

Note

p is not accessed.

solveFwdIncremental(sol, rhs)[source]

Solve the linearized (incremental) forward problem for a given right-hand side

Parameters:
  • sol (dolfin.Vector) – Solution of the incremental forward problem (state)

  • rhs (dolfin.Vector) – Right hand side of the linear system

soupy.modeling.controlQoI

class soupy.modeling.controlQoI.ControlQoI[source]

Bases: object

Abstract class to define the optimization quantity of interest for the optimal control problem under uncertainty In the following x will denote the variable [u, m, p, z], denoting respectively the state u, the parameter m, the adjoint variable p, and the control variable z

The methods in the class ControlQoI will usually access the state u and possibly the parameter m and control z. The adjoint variables will never be accessed.

apply_ij(i, j, dir, out)[source]

Apply the second variation \(\delta_{ij}\) (i,j = soupy.STATE, soupy.PARAMETER, soupy.CONTROL) of the cost in direction dir.

cost(x)[source]

Given x evaluate the cost functional. Only the state u and (possibly) the parameter m are accessed.

grad(i, x, out)[source]

Given the state and the paramter in x, compute the partial gradient of the QoI in with respect to the state (i == soupy.STATE), parameter (i == soupy.PARAMETER), or control (i == soupy.CONTROL).

setLinearizationPoint(x, gauss_newton_approx=False)[source]

Set the point for linearization.

class soupy.modeling.controlQoI.L2MisfitControlQoI(Vh, ud)[source]

Bases: ControlQoI

Class for the \(L^2(\Omega)\) misfit functional,

\[\int_\Omega (u - u_d)^2 dx\]

where \(u_d\) is the reference state.

Constructor

Parameters:
  • Vh (list of dolfin.FunctionSpace) – List of function spaces for the state, parameter, adjoint, and optimization variables

  • ud (dolfin.Vector) – The reference state as a vector

adj_rhs(x, rhs)[source]

The right hand for the adjoint problem (i.e. the derivative of the Lagrangian functional with respect to the state u).

Parameters:
  • x (list of dolfin.Vector) – List of vectors [u, m, p, z] representing the state, parameter, adjoint, and control variables

  • rhs (dolfin.Vector) – The assembled rhs for the adjoint problem.

apply_ij(i, j, dir, out)[source]

Apply the second variation \(\delta_{ij}\) (i,j = soupy.STATE, soupy.PARAMETER, soupy.CONTROL) of the QoI in direction dir.

Parameters:
  • i (int) – Index of the output variable

  • j (int) – Index of the input variable

  • dir (dolfin.Vector) – The direction in which to apply the second variation

  • out (dolfin.Vector) – The assembled second variation

..note:: setLinearizationPoint must be called before calling this method.

apply_ijk(i, j, k, dir1, dir2, out)[source]

Apply the third order variation of the QoI in the i th, j th, and k th variables in directions dir1 and dir2.

Parameters:
  • i (int) – First variable index (soupy.STATE, soupy.PARAMETER, or soupy.CONTROL)

  • j (int) – Second variable index (soupy.STATE, soupy.PARAMETER, or soupy.CONTROL)

  • k (int) – Third variable index (soupy.STATE, soupy.PARAMETER, or soupy.CONTROL)

  • dir1 (dolfin.Vector) – Direction for variable j

  • dir2 (dolfin.Vector) – Direction for variable k

  • out (dolfin.Vector) – The assembled third variation

cost(x)[source]

Evaluate the qoi at given point \(q(u,m,z)\)

Parameters:

x (list of dolfin.Vector) – List of vectors [u, m, p, z] representing the state, parameter, adjoint, and control variables

Returns:

QoI evaluated at x

grad(i, x, out)[source]

First variation of the QoI with respect to the i th variable where i is either soupy.STATE, soupy.PARAMETER, or soupy.CONTROL.

Parameters:
  • i (int) – Index of the variable with respect to which the first variation is taken

  • x (list of dolfin.Vector) – List of vectors [u, m, p, z] representing the state, parameter, adjoint, and control variables

  • out (dolfin.Vector) – The assembled first variation

setLinearizationPoint(x, gauss_newton_approx=False)[source]

Specify the linearization point for computation of the second variations in method apply_ij.

Parameters:

x (list of dolfin.Vector) – List of vectors [u, m, p, z] representing the state, parameter, adjoint, and control variables

class soupy.modeling.controlQoI.L2MisfitVarfHandler(ud, chi=None)[source]

Bases: object

Form handler for the \(L^2\) Misfit

\[\int_{\Omega} \chi (u - u_d)^2 dx\]

where \(u_d\) is the reference state and \(\chi\) is the characteristic function defining the region of integration

Constructor

Parameters:
  • ud (dolfin.Function or dolfin.Expression) – The reference state

  • chi (dolfin.Function or dolfin.Expression) – The characteristic function defining the region of integration

class soupy.modeling.controlQoI.VariationalControlQoI(Vh, form_handler)[source]

Bases: ControlQoI

Class for a QoI defined by its variational form

Constructor

Parameters:
  • Vh (list of dolfin.FunctionSpace) – List of function spaces for the state, parameter, adjoint, and optimization variables

  • form_handler – The form handler for the variational form with a __call__ method that takes as input the state, parameter, and control variables as functions and returns the variational form

adj_rhs(x, rhs)[source]

The right hand for the adjoint problem (i.e. the derivative of the Lagrangian functional with respect to the state u).

Parameters:
  • x (list of dolfin.Vector) – List of vectors [u, m, p, z] representing the state, parameter, adjoint, and control variables

  • rhs (dolfin.Vector) – The assembled rhs for the adjoint problem.

apply_ij(i, j, dir, out)[source]

Apply the second variation \(\delta_{ij}\) (i,j = soupy.STATE, soupy.PARAMETER, soupy.CONTROL) of the QoI in direction dir.

Parameters:
  • i (int) – Index of the output variable

  • j (int) – Index of the input variable

  • dir (dolfin.Vector) – The direction in which to apply the second variation

  • out (dolfin.Vector) – The assembled second variation

..note:: setLinearizationPoint must be called before calling this method.

apply_ijk(i, j, k, dir1, dir2, out)[source]

Apply the third order variation of the QoI in the i th, j th, and k th variables in directions dir1 and dir2.

Parameters:
  • i (int) – First variable index (soupy.STATE, soupy.PARAMETER, or soupy.CONTROL)

  • j (int) – Second variable index (soupy.STATE, soupy.PARAMETER, or soupy.CONTROL)

  • k (int) – Third variable index (soupy.STATE, soupy.PARAMETER, or soupy.CONTROL)

  • dir1 (dolfin.Vector) – Direction for variable j

  • dir2 (dolfin.Vector) – Direction for variable k

  • out (dolfin.Vector) – The assembled third variation

cost(x)[source]

Evaluate the qoi at given point \(q(u,m,z)\)

Parameters:

x (list of dolfin.Vector) – List of vectors [u, m, p, z] representing the state, parameter, adjoint, and control variables

Returns:

QoI evaluated at x

grad(i, x, out)[source]

First variation of the QoI with respect to the i th variable where i is either soupy.STATE, soupy.PARAMETER, or soupy.CONTROL.

Parameters:
  • i (int) – Index of the variable with respect to which the first variation is taken

  • x (list of dolfin.Vector) – List of vectors [u, m, p, z] representing the state, parameter, adjoint, and control variables

  • out (dolfin.Vector) – The assembled first variation

setLinearizationPoint(x, gauss_newton_approx=False)[source]

Specify the linearization point for computation of the second variations in method apply_ij.

Parameters:

x (list of dolfin.Vector) – List of vectors [u, m, p, z] representing the state, parameter, adjoint, and control variables

soupy.modeling.controlModelHessian

class soupy.modeling.controlModelHessian.ControlModelHessian(model)[source]

Bases: object

Hessian of the the QoI map with respect to the CONTROL variable \(z\) wrapping around a soupy.ControlModel

Constructor:

Parameters:

model (soupy.ControlModel) – Control model defining the control to QoI map

init_vector(z, dim)[source]

Initialize vector to be compatible with operator. Note since the Hessian is symmetric, dimension is not important

mult(zhat, Hzhat)[source]

Apply the Hessian of the QoI.

Parameters:
  • zhat (dolfin.Vector or similar) – Direction for Hessian action

  • Hzhat (dolfin.Vector or similar) – Output for Hessian action

::note Assumes the model.setLinearizationPoint has been called

soupy.modeling.controlCostFunctional

class soupy.modeling.controlCostFunctional.ControlCostFunctional[source]

Bases: object

Base class for the cost function for solving an optimal control problem under uncertainty.

cost(z, order=0)[source]

Given the control variable z evaluate the cost functional. Order specifies the order of derivatives required after the computation of the cost

generate_vector(component='ALL')[source]

If component is STATE, PARAMETER, ADJOINT, or CONTROL, return a vector corresponding to that function space. If component is "ALL", Generate the list of vectors x = [u,m,p,z]

grad(g)[source]

Evaluate the gradient of the cost functional. Assumes cost is called with order >=1

hessian(zhat, Hzhat)[source]

Evaluate the Hessian of the cost functional acting in direction zhat. Assumes cost is called with order >=2

class soupy.modeling.controlCostFunctional.DeterministicControlCostFunctional(model, prior, penalization=None)[source]

Bases: ControlCostFunctional

This class implements a deterministic approximation for the optimal control problem under uncertainty by considering a fixed parameter at the mean of the prior

\[J(z) := Q(\bar{m}, z) + P(z)\]

Constructor

Parameters:
  • model (soupy.ControlModel) – control model containing the soupy.PDEVariationalControlProblem and soupy.ControlQoI

  • prior (hippylib.Prior) – The prior distribution for the random parameter

  • penalization (soupy.Penalization) – An optional penalization object for the cost of control

computeComponents(z, order=0)[source]

Computes the components for the evaluation of the cost functional

Parameters:
  • z (dolfin.Vector) – the control variable

  • order (int) – Order of the derivatives needed. 0 for cost, 1 for gradient, 2 for Hessian

cost(z, order=0, **kwargs)[source]

Computes the value of the cost functional at given control \(z\)

Parameters:
  • z (dolfin.Vector) – the control variable

  • order (int) – Order of the derivatives needed after evaluation 0 for cost, 1 for gradient, 2 for Hessian

Returns:

Value of the cost functional

generate_vector(component='ALL')[source]

If component is STATE, PARAMETER, ADJOINT, or CONTROL, return a vector corresponding to that function space. If component is "ALL", Generate the list of vectors x = [u,m,p,z]

grad(g)[source]

Computes the gradient of the cost functional

Parameters:

g (dolfin.Vector) – (Dual of) the gradient of the cost functional

Returns:

the norm of the gradient in the correct inner product \((g_z,g_z)_{Z}^{1/2}\)

Note

Assumes self.cost has been called with order >= 1

hessian(zhat, Hzhat)[source]

Apply the the reduced Hessian to the vector \(zhat\)

Parameters:
  • zhat (dolfin.Vector) – The direction of Hessian action

  • Hzhat (dolfin.Vector) – The assembled Hessian action

Note

Assumes self.cost has been called with order >= 2

objective(z)[source]
class soupy.modeling.controlCostFunctional.PenalizationControlCostFunctional(Vh, penalization, augment_control=False)[source]

Bases: ControlCostFunctional

Cost functional that is defined by a soupy.Penalization object. Used for optimization problems where the cost depends only on the control variable \(z\) and not the PDE.

Constructor:

Parameters:
  • Vh – List of function spaces for the state, parameter, adjoint, and optimization variables

  • penalization – Penalization term defining the cost functional

  • augment_control (bool) – Set to True if optimization problem uses a soupy.AugmentedVector as the control/optimization variable.

cost(z, order=0)[source]

Given the control variable z evaluate the cost functional. Order specifies the order of derivatives required after the computation of the cost

generate_vector(component='ALL')[source]

If component is STATE, PARAMETER, ADJOINT, CONTROL return a vector corresponding to that function space. If component == CONTROL and self.augment_control is True, will return an soupy.AugmentedVector that augments the control variable z with a scalar that can be used for optimization. If component == "ALL", Generate the list of vectors x = [u,m,p,z]. Note that in this case, the CONTROL variable will not be augmented with the scalar, and can be used directly for methods like solveFwd.

grad(g)[source]

Evaluate the gradient of the cost functional. Assumes cost is called with order >=1

hessian(zhat, Hzhat)[source]

Evaluate the Hessian of the cost functional acting in direction zhat. Assumes cost is called with order >=2

class soupy.modeling.controlCostFunctional.RiskMeasureControlCostFunctional(risk_measure, penalization=None)[source]

Bases: object

This class implements a risk measure cost functional for optimal control problem under uncertainty

\[J(z) := \rho[Q(m, z)] + P(z)\]

Constructor

Parameters:
  • risk_measure (soupy.RiskMeasure) – Class implementing the risk measure \(\rho(m,z)\)

  • penalization (soupy.Penalization) – An optional penalization object for the cost of control

cost(z, order=0, **kwargs)[source]

Computes the value of the cost functional at given control \(z\)

Parameters:
  • z (dolfin.Vector) – the control variable

  • order (int) – Order of the derivatives needed after evaluation 0 for cost, 1 for gradient, 2 for Hessian

  • kwargs – additional arguments, e.g. rng for the risk measure computation

Returns:

Value of the cost functional

generate_vector(component='ALL')[source]

If component is STATE, PARAMETER, ADJOINT, CONTROL return a vector corresponding to that function space. If component == CONTROL and the underlying risk measure uses a soupy.AugmentedVector, return an soupy.AugmentedVector that augments the control variable z with a scalar that can be used for optimization. If component == "ALL", Generate the list of vectors x = [u,m,p,z]. Note that in this case, the CONTROL variable will not be augmented with the scalar, and can be used directly for methods like solveFwd.

grad(g)[source]

Computes the gradient of the cost functional

Parameters:

g (dolfin.Vector) – (Dual of) the gradient of the cost functional

Returns:

the norm of the gradient in the correct inner product \((g_z,g_z)_{Z}^{1/2}\)

Note

Assumes self.cost has been called with order >= 1

hessian(zhat, Hzhat)[source]

Apply the the reduced Hessian to the vector \(zhat\)

Parameters:
  • zhat (dolfin.Vector) – The direction of Hessian action

  • Hzhat (dolfin.Vector) – The assembled Hessian action

Note

Assumes self.cost has been called with order >= 2

soupy.modeling.controlCostHessian

class soupy.modeling.controlCostHessian.ControlCostHessian(cost_functional)[source]

Bases: object

Hessian of the the QoI map with respect to the CONTROL variable \(z\) wrapping around a soupy.ControlModel. Hessian will be applied where soupy.ControlModel is last called with order>=2

Constructor:

Parameters:

cost_functional (soupy.ControlCostFunctional) – cost functional

inner(z1, z2)[source]

Computes the Hessian weighted inner product

..math:: z_1^T H z_2

Parameters:
  • z1 (dolfin.Vector or similar) – Vector 1

  • z2 (dolfin.Vector or similar) – Vector 2

Returns:

The Hessian weighted inner product

mult(zhat, Hzhat)[source]

Apply the Hessian of the QoI in direction \(\hat{z}\)

Parameters:
  • zhat (dolfin.Vector or similar) – Direction for Hessian action

  • Hzhat (dolfin.Vector or similar) – Output for Hessian action

property ncalls

soupy.modeling.riskMeasure

class soupy.modeling.riskMeasure.RiskMeasure[source]

Bases: object

Abstract class for the risk measure \(\rho[Q](z)\)

computeComponents(z, order=0)[source]

Computes the components for the evaluation of the risk measure

Parameters:
  • z (dolfin.Vector) – the control variable

  • order (int) – Order of the derivatives needed. 0 for cost, 1 for gradient, 2 for Hessian

cost()[source]

Evaluates the value of the risk measure once components have been computed

Returns:

Value of the cost functional

Note

Assumes computeComponents has been called with order>=0

generate_vector(components='ALL')[source]

If components is STATE, PARAMETER, ADJOINT, or CONTROL, return a vector corresponding to that function space. If components is "ALL", Generate the list of vectors x = [u,m,p,z]

grad()[source]

Evaluates the gradient of the risk measure once components have been computed

Parameters:

g (dolfin.Vector) – (Dual of) the gradient of the risk measure to store result in

Note

Assumes self.computeComponents has been called with order >= 1

hessian(zhat, Hzhat)[source]

Apply the hessian of the risk measure once components have been computed in the direction zhat

Parameters:
  • zhat (dolfin.Vector) – Direction for application of Hessian action of the risk measure

  • Hzhat (dolfin.Vector) – (Dual of) Result of the Hessian action of the risk measure to store the result in

Note

Assumes self.computeComponents has been called with order >= 1

soupy.modeling.meanVarRiskMeasureStochastic

class soupy.modeling.meanVarRiskMeasureStochastic.MeanVarRiskMeasureStochastic(control_model, prior, settings=hippylib.ParameterList)[source]

Bases: RiskMeasure

Class for memory efficient evaluation of the Mean + Variance risk measure

\[\rho[Q](z) = \mathbb{E}_m[Q(m,z)] + \beta \mathbb{V}_m[Q(m,z)]\]

that allows for stochastic approximations and SGD

Constructor:

Parameters:
  • control_model (soupy.ControlModel) – control model containing the soupy.PDEVariationalControlProblem and soupy.ControlQoI

  • prior (hippylib.Prior) – The prior distribution for the random parameter

  • settings (hippylib.ParameterList) – additional settings

computeComponents(z, order=0, sample_size=100, rng=None)[source]

Computes the components for the evaluation of the risk measure

Parameters:
  • z (dolfin.Vector) – the control variable

  • order (int) – Order of the derivatives needed. 0 for cost, 1 for gradient, 2 for Hessian

  • rng (hippylib.Random) – The random number generator used for sampling. If None the default uses hippylib.parRandom

cost()[source]

Evaluates the value of the risk measure once components have been computed

Returns:

Value of the cost functional

Note

Assumes computeComponents has been called with order>=0

generate_vector(component='ALL')[source]

If components is STATE, PARAMETER, ADJOINT, or CONTROL, return a vector corresponding to that function space. If components is "ALL", Generate the list of vectors x = [u,m,p,z]

grad(g)[source]

Evaluates the gradient of the risk measure once components have been computed

Parameters:

g (dolfin.Vector) – (Dual of) the gradient of the risk measure to store result in

Note

Assumes self.computeComponents has been called with order >= 1

hessian(zhat, Hzhat)[source]

Apply the hessian of the risk measure once components have been computed in the direction zhat

Parameters:
  • zhat (dolfin.Vector) – Direction for application of Hessian action of the risk measure

  • Hzhat (dolfin.Vector) – (Dual of) Result of the Hessian action of the risk measure to store the result in

Note

Assumes self.computeComponents has been called with order >= 1

soupy.modeling.meanVarRiskMeasureStochastic.meanVarRiskMeasureStochasticSettings(data={'beta': [0, 'Weighting factor for variance']})[source]

soupy.modeling.meanVarRiskMeasureSAA

class soupy.modeling.meanVarRiskMeasureSAA.MeanVarRiskMeasureSAA(control_model, prior, settings=hippylib.ParameterList, comm_sampler=mpi4py.MPI.COMM_WORLD)[source]

Bases: RiskMeasure

Mean + variance risk measure using sample average approximation

\[\rho[Q](z) = \mathbb{E}_m[Q(m,z)] + \beta \mathbb{V}_m[Q(m,z)]\]

with sample parallelism using MPI

Note

currently does not support simultaneous sample and mesh partition parallelism

Constructor:

Parameters:
  • control_model (soupy.ControlModel) – control model containing the soupy.PDEVariationalControlProblem and soupy.ControlQoI

  • prior (hippylib.Prior) – The prior distribution for the random parameter

  • settings (hippylib.ParameterList) – additional settings

  • comm_sampler (mpi4py.MPI.Comm) – MPI communicator for sample parallelism

computeComponents(z, order=0, **kwargs)[source]

Computes the components for the evaluation of the risk measure

Parameters:
  • z (dolfin.Vector) – the control variable

  • order (int) – Order of the derivatives needed. 0 for cost, 1 for gradient, 2 for Hessian

  • kwargs – dummy keyword arguments for compatibility

cost()[source]

Evaluates the value of the risk measure once components have been computed

Returns:

Value of the cost functional

Note

Assumes computeComponents has been called with order>=0

gather_samples()[source]

Gather the QoI samples from all processes

Returns:

An array of the sample QoI values

Return type:

numpy.ndarray

generate_vector(component='ALL')[source]

If components is STATE, PARAMETER, ADJOINT, or CONTROL, return a vector corresponding to that function space. If components is "ALL", generate the list of vectors x = [u,m,p,z]

grad(g)[source]

Evaluates the gradient of the risk measure once components have been computed

Parameters:

g (dolfin.Vector) – (Dual of) the gradient of the risk measure to store result in

Note

Assumes self.computeComponents has been called with order >= 1

hessian(zhat, Hzhat)[source]

Apply the hessian of the risk measure once components have been computed in the direction zhat

Parameters:
  • zhat (dolfin.Vector) – Direction for application of Hessian action of the risk measure

  • Hzhat (dolfin.Vector) – (Dual of) Result of the Hessian action of the risk measure to store the result in

Note

Assumes self.computeComponents has been called with order >= 1

soupy.modeling.meanVarRiskMeasureSAA.meanVarRiskMeasureSAASettings(data={'beta': [0, 'Weighting factor for variance'], 'sample_size': [100, 'Number of Monte Carlo samples'], 'seed': [1, 'rng seed for sampling']})[source]

soupy.modeling.augmentedVector

class soupy.modeling.augmentedVector.AugmentedVector(v, copy_vector=True)[source]

Bases: object

Class representing an augmented optimization variable \((z, t)\), where \(t\) is a real number. Methods mirror that of dolfin.Vector. Assumes the last element in the array is the t variable

Note

Currently supports only sample parallelism (i.e. serial mesh)

Parameters:
  • vdolfin.Vector to be augmented

  • copy_vector – If True, copy the vector, otherwise use the same memory

add_local(vt_array)[source]
apply(method)[source]
axpy(a, vt)[source]
copy()[source]
get_local()[source]
get_scalar()[source]
get_vector()[source]
inner(vt)[source]
set_local(vt_array)[source]
set_scalar(t)[source]
set_vector(v)[source]
zero()[source]

soupy.modeling.superquantileRiskMeasureSAA

class soupy.modeling.superquantileRiskMeasureSAA.SuperquantileRiskMeasureSAA(control_model, prior, settings=hippylib.ParameterList, comm_sampler=mpi4py.MPI.COMM_WORLD)[source]

Bases: RiskMeasure

Risk measure for the sample average approximation of the superquantile risk measure (CVaR) with sample parallelism using MPI

Constructor:

Parameters:
  • control_model (soupy.ControlModel) – control model containing the soupy.PDEVariationalControlProblem and soupy.ControlQoI

  • prior (hippylib.Prior) – The prior distribution for the random parameter

  • settings (hippylib.ParameterList) – additional settings given in superquantileRiskMeasureSAASettings()

  • comm_sampler (mpi4py.MPI.Comm) – MPI communicator for sample parallelism

computeComponents(zt, order=0, **kwargs)[source]

Computes the components for the evaluation of the risk measure

Parameters:
  • zt (soupy.AugmentedVector) – the control variable with the scalar t appended

  • order (int) – Order of the derivatives needed. 0 for cost, 1 for gradient, 2 for Hessian

  • kwargs – dummy keyword arguments for compatibility

cost()[source]

Evaluates the value of the risk measure once components have been computed

Returns:

Value of the cost functional

Note

Assumes computeComponents has been called with order>=0

gather_samples()[source]

Gather the QoI samples from all processes

Returns:

An array of the sample QoI values

Return type:

numpy.ndarray

generate_vector(component='ALL')[source]

If component is STATE, PARAMETER, ADJOINT, return a vector corresponding to that function space. If component == CONTROL, return an soupy.AugmentedVector that augments the control variable z with a scalar that can be used for optimization

If component == "ALL", Generate the list of vectors x = [u,m,p,z]. Note that in this case, the CONTROL variable will not be augmented with the scalar, and can be used directly for methods like solveFwd.

grad(gt)[source]

Evaluates the gradient of the risk measure once components have been computed

Parameters:

g (dolfin.Vector) – (Dual of) the gradient of the risk measure to store result in

Note

Assumes self.computeComponents has been called with order>=1

hessian(zt_hat, Hzt_hat)[source]

Apply the hessian of the risk measure once components have been computed in the direction zhat

Parameters:
  • zt_hat (soupy.AugmentedVector) – Direction for application of Hessian action of the risk measure

  • Hzt_hat (soupy.AugmentedVector) – (Dual of) Result of the Hessian action of the risk measure to store the result in

Note

Assumes self.computeComponents has been called with order >= 1

superquantile()[source]

Evaluate the superquantile using the computed samples

Returns:

Value of the superquantile by sampling

Note

Assumes computeComponents has been called with order>=0

soupy.modeling.superquantileRiskMeasureSAA.sample_superquantile(samples, beta)[source]

Evaluate superquantile from samples

soupy.modeling.superquantileRiskMeasureSAA.sample_superquantile_by_minimization(samples, beta, epsilon=0.01)[source]

Evaluate superquantile from samples by minimization

soupy.modeling.superquantileRiskMeasureSAA.superquantileRiskMeasureSAASettings(data={'beta': [0.95, 'Quantile value for superquantile'], 'epsilon': [0.01, 'Sharpness of smooth plus approximation'], 'sample_size': [100, 'Number of Monte Carlo samples'], 'seed': [1, 'rng seed for sampling'], 'smoothplus_type': ['quartic', 'approximation type for smooth plus function']})[source]

soupy.modeling.transformedMeanRiskMeasureSAA

class soupy.modeling.transformedMeanRiskMeasureSAA.FunctionWrapper(function, grad=None, hessian=None, delta=0.0001)[source]

Bases: object

Wrapper for a scalar function \(f(x)\) that can be called with a single argument. Will use finite difference for gradient and Hessian if not provided.

Constructor:

Parameters:
  • function – Callable function for \(f(x)\)

  • grad – Optional callable function for the gradient

  • hessian – Optional callable function for the Hessian

  • delta – Finite difference step for derivatives if not provided

grad(x)[source]
hessian(x)[source]
class soupy.modeling.transformedMeanRiskMeasureSAA.IdentityFunction[source]

Bases: object

Scalar function \(f(x) = x\). Supports both x as a single scalar or a vector of scalars.

grad(x)[source]
hessian(x)[source]
class soupy.modeling.transformedMeanRiskMeasureSAA.TransformedMeanRiskMeasureSAA(control_model, prior, settings=hippylib.ParameterList, comm_sampler=mpi4py.MPI.COMM_WORLD)[source]

Bases: RiskMeasure

Risk measure defined by transformations of an expectation of the following form

\[\rho[Q](z) = g \left( \mathbb{E}_m[f(Q(m,z))] \right)\]

given user-defined inner/outer scalar functions \(f(x)\) and \(g(x)\), with sample parallelism using MPI.

Note

currently does not support simultaneous sample and mesh partition parallelism

Constructor:

Parameters:
  • control_model (soupy.ControlModel) – control model containing the soupy.PDEVariationalControlProblem and soupy.ControlQoI

  • prior (hippylib.Prior) – The prior distribution for the random parameter

  • settings (hippylib.ParameterList) – Additional settings

  • comm_sampler (mpi4py.MPI.Comm) – MPI communicator for sample parallelism

computeComponents(z, order=0, **kwargs)[source]

Computes the components for the evaluation of the risk measure

Parameters:
  • z (dolfin.Vector) – the control variable

  • order (int) – Order of the derivatives needed. 0 for cost, 1 for gradient, 2 for Hessian

  • kwargs – dummy keyword arguments for compatibility

cost()[source]

Evaluates the value of the risk measure once components have been computed

Returns:

Value of the cost functional

Note

Assumes computeComponents has been called with order>=0

gather_samples()[source]

Gather the QoI samples from all processes

Returns:

An array of the sample QoI values

Return type:

numpy.ndarray

generate_vector(component='ALL')[source]

If component is STATE, PARAMETER, ADJOINT, or CONTROL, return a vector corresponding to that function space. If component is "ALL", Generate the list of vectors x = [u,m,p,z]

grad(g)[source]

Evaluates the gradient of the risk measure once components have been computed

Parameters:

g (dolfin.Vector) – (Dual of) the gradient of the risk measure to store result in

Note

Assumes self.computeComponents has been called with order >= 1

hessian(zhat, Hzhat)[source]

Apply the hessian of the risk measure once components have been computed in the direction zhat

Parameters:
  • zhat (dolfin.Vector) – Direction for application of Hessian action of the risk measure

  • Hzhat (dolfin.Vector) – (Dual of) Result of the Hessian action of the risk measure to store the result in

Note

Assumes self.computeComponents has been called with order >= 1

soupy.modeling.transformedMeanRiskMeasureSAA.transformedMeanRiskMeasureSAASettings(data={'inner_function': [<soupy.modeling.transformedMeanRiskMeasureSAA.IdentityFunction object>, 'Transformation of Q inside the expectation'], 'outer_function': [<soupy.modeling.transformedMeanRiskMeasureSAA.IdentityFunction object>, 'Transformation of Q outside the expectation'], 'sample_size': [100, 'Number of Monte Carlo samples'], 'seed': [1, 'rng seed for sampling']})[source]

soupy.modeling.penalization

class soupy.modeling.penalization.L2Penalization(Vh, alpha)[source]

Bases: Penalization

\(L^2(\Omega)\) integral over the domain

\[P(z) = \alpha \int_{\Omega} |z|^2 dx\]

For finite dimensional controls z, this amounts to a little \(\ell_2\) norm In this case, Vh[soupy.CONTROL] needs to be a dolfin.VectorFunctionSpace of reals

Constructor

Parameters:
  • Vh (list of dolfin.FunctionSpace) – List of function spaces the state, parameter, adjoint, and control

  • alpha (float) – Weight for the penalization term

cost(z)[source]

Compute the penalization at given control \(z\)

Parameters:

z (dolfin.Vector or AugmentedVector) – The control variable

grad(z, out)[source]

Compute the gradient of the penalty at given control \(z\)

Parameters:
  • z (dolfin.Vector or AugmentedVector) – The control variable

  • out (dolfin.Vector or AugmentedVector) – The assembled gradient vector

hessian(z, zhat, out)[source]

Compute the Hessian of the penalty at given control \(z\)

Parameters:
  • z (dolfin.Vector or soupy.AugmentedVector) – The control variable

  • zhat (dolfin.Vector or soupy.AugmentedVector) – The direction for Hessian action

  • out (dolfin.Vector or soupy.AugmentedVector) – The assembled Hessian action vector

init_vector(v, dim=0)[source]
class soupy.modeling.penalization.MultiPenalization(Vh, penalization_list, alpha_list=None)[source]

Bases: Penalization

Penalization term for the sum of individual penalties

\[P(z) = \sum_{i=1}^{n} \alpha_i P_i(z)\]

Constructor

Parameters:
  • Vh (list of dolfin.FunctionSpace) – List of function spaces the state, parameter, adjoint, and control

  • penalization_list (list of Penalization) – List of penalization objects

  • alpha_list (list of floats) – List of weights for each penalization term

cost(z)[source]

Compute the penalization at given control \(z\)

Parameters:

z (dolfin.Vector or soupy.AugmentedVector) – The control variable

grad(z, out)[source]

Compute the gradient of the penalty at given control \(z\)

Parameters:
  • z (dolfin.Vector or soupy.AugmentedVector) – The control variable

  • out (dolfin.Vector or soupy.AugmentedVector) – The assembled gradient vector

hessian(z, zhat, out)[source]

Compute the Hessian of the penalty at given control \(z\)

Parameters:
  • z (dolfin.Vector or soupy.AugmentedVector) – The control variable

  • zhat (dolfin.Vector or soupy.AugmentedVector) – The direction for Hessian action

  • out (dolfin.Vector or soupy.AugmentedVector) – The assembled Hessian action vector

init_vector(z)[source]
class soupy.modeling.penalization.Penalization[source]

Bases: object

Abstract class for the penalization on the control variable \(z\)

cost(z)[source]
grad(z, out)[source]
hessian(z, zhat, out)[source]
init_vector(z)[source]
class soupy.modeling.penalization.VariationalPenalization(Vh, form_handler)[source]

Bases: Penalization

Penalization given by a variational form in terms of the control variable \(z\)

Constructor:

Parameters:
  • Vh (list of dolfin.FunctionSpace) – List of function spaces for the state, parameter, adjoint, and optimization variables

  • form_handler – The form handler for the penalization with a __call__ method that takes the control variable z as a dolfin.Function and returns the variational form.

cost(z)[source]
grad(z, out)[source]
hessian(z, zhat, out)[source]
init_vector(z)[source]
class soupy.modeling.penalization.WeightedL2Penalization(Vh, M, alpha)[source]

Bases: Penalization

A weighted L2 norm penalization
\[P(z) = z^T M z\]

where \(M\) is a symmetric positive definite weight matrix

Constructor:

Parameters:
  • Vh (list of dolfin.FunctionSpace) – List of function spaces the state, parameter, adjoint, and control

  • M (dolfin.Matrix like) – Weighting matrix with method mult

  • alpha (float) – Weight for the penalization term

cost(z)[source]

Compute the penalization at given control \(z\)

Parameters:

z (dolfin.Vector or soupy.AugmentedVector) – The control variable

grad(z, out)[source]

Compute the gradient of the penalty at given control \(z\)

Parameters:
  • z (dolfin.Vector or soupy.AugmentedVector) – The control variable

  • out (dolfin.Vector or soupy.AugmentedVector) – The assembled gradient vector

hessian(z, zhat, out)[source]

Compute the Hessian of the penalty at given control \(z\)

Parameters:
  • z (dolfin.Vector or soupy.AugmentedVector) – The control variable

  • zhat (dolfin.Vector or soupy.AugmentedVector) – The direction for Hessian action

  • out (dolfin.Vector or soupy.AugmentedVector) – The assembled Hessian action vector

init_vector(v, dim=0)[source]

soupy.modeling.variables

Enumerator for the variables of the inverse problem:
  • STATE = 0 for the state variable

  • PARAMETER = 1 for the parameter variable

  • ADJOINT = 2 for the adjoint variable

  • CONTROL = 3 for the control variable