soupy.modeling
soupy.modeling.PDEControlProblem
- class soupy.modeling.PDEControlProblem.PDEVariationalControlProblem(*args: Any, **kwargs: Any)[source]
Bases:
PDEVariationalProblemConstructor
- Parameters:
Vh (list of
dolfin.FunctionSpace) – List of function spaces the state, parameter, adjoint, and controlvarf_handler – Variational form handler with
__call__methodbc – 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}\).
- 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}.\)
- 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:
objectThis class contains the structure needed to evaluate the control objective As inputs it takes a
PDEVariationalControlProblem, and aQoiobject.- In the following we will denote with
uthe state variablemthe (model) parameter variablepthe adjoint variablezthe control variable
Constructor
- Parameters:
problem (
PDEVariationalControlProblem) – The PDE problemqoi (
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 variableout (
dolfin.Vector) – The action of the \(C_z\) block ondm
Note
This routine assumes that
outhas 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 variableout (
dolfin.Vector) – The action of the \(C_{m}^T\) block ondp
..note:: This routine assumes that
outhas 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 variableout (
dolfin.Vector) – The action of the \(C_z\) block ondz
Note
This routine assumes that
outhas 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 variableout (
dolfin.Vector) – The action of the \(C_{z}^T\) block ondp
..note:: This routine assumes that
outhas 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 variableout (
dolfin.Vector) – The action of the \(W_{mm}\) block ondu
Note
This routine assumes that
outhas 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 variableout (
dolfin.Vector) – The action of the \(W_{mu}\) block ondu
Note
This routine assumes that
outhas 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 variableout (
dolfin.Vector) – The action of the \(W_{um}\) block ondu
Note
This routine assumes that
outhas 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 variableout (
dolfin.Vector) – The action of the \(W_{uu}\) block ondu
Note
This routine assumes that
outhas 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 variableout (
dolfin.Vector) – The action of the \(W_{uz}\) block ondu
Note
This routine assumes that
outhas 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 variableout (
dolfin.Vector) – The action of the \(W_{zu}\) block ondu
Note
This routine assumes that
outhas 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 variableout (
dolfin.Vector) – The action of the \(W_{zz}\) block ondu
Note
This routine assumes that
outhas 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, orsoupy.CONTROLj (int) – The input variable index
soupy.STATE,soupy.PARAMETER,soupy.ADJOINT, orsoupy.CONTROLd (
dolfin.Vector) – The vector to which the Hessian is appliedout (
dolfin.Vector) – The action of the Hessian ond
- cost(x)[source]
Evaluate the QoI at the a given point \(q(u,m,z)\).
- Parameters:
x (list of
dolfin.Vectorobjects) – The pointx = [u,m,p,z]at which to evaluate the QoI- Returns:
The value of the QoI
Note
pis 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.Vectorobjects) – The pointx = [u,m,p,z]. Provides the state, parameter, adjoint, and control variablesu,m,p, andzfor the assembling the gradient action Vectoruis 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.Vectorobjects) – The pointx = [u,m,p,z]. Provides the state, parameter, adjoint, and control variablesu,m,p, andzfor the assembling the gradient action. Vectoruis 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 adolfin.Vectorobject of the appropriate size Ifcomponent = soupy.STATEreturn onlyu. Ifcomponent = soupy.PARAMETERreturn onlym. Ifcomponent = soupy.ADJOINTreturn onlyp. Ifcomponent = soupy.CONTROLreturn onlyz.
- 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 evaluatedgauss_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 variablesu,m, andzfor the solution assembling the adjoint operator. Vectoruis also used to assemble the adjoint right hand side.
- Type:
list of
dolfin.Vectorobjects
Note
pis 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 variablesmandzfor the solution of the forward problem and the initial guessuif the forward problem is non-linear
- Type:
list of
dolfin.Vectorobjects
Note
pis not accessed.
soupy.modeling.controlQoI
- class soupy.modeling.controlQoI.ControlQoI[source]
Bases:
objectAbstract class to define the optimization quantity of interest for the optimal control problem under uncertainty In the following
xwill denote the variable[u, m, p, z], denoting respectively the stateu, the parameterm, the adjoint variablep, and the control variablezThe methods in the class ControlQoI will usually access the state
uand possibly the parametermand controlz. 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 directiondir.
- cost(x)[source]
Given
xevaluate the cost functional. Only the stateuand (possibly) the parametermare accessed.
- class soupy.modeling.controlQoI.L2MisfitControlQoI(Vh, ud)[source]
Bases:
ControlQoIClass 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 variablesud (
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 variablesrhs (
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 directiondir.- 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 variationout (
dolfin.Vector) – The assembled second variation
..note::
setLinearizationPointmust 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
ith,jth, andkth variables in directionsdir1anddir2.- Parameters:
i (int) – First variable index (
soupy.STATE,soupy.PARAMETER, orsoupy.CONTROL)j (int) – Second variable index (
soupy.STATE,soupy.PARAMETER, orsoupy.CONTROL)k (int) – Third variable index (
soupy.STATE,soupy.PARAMETER, orsoupy.CONTROL)dir1 (
dolfin.Vector) – Direction for variablejdir2 (
dolfin.Vector) – Direction for variablekout (
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
ith variable whereiis eithersoupy.STATE,soupy.PARAMETER, orsoupy.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 variablesout (
dolfin.Vector) – The assembled first variation
- class soupy.modeling.controlQoI.L2MisfitVarfHandler(ud, chi=None)[source]
Bases:
objectForm 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.Functionordolfin.Expression) – The reference statechi (
dolfin.Functionordolfin.Expression) – The characteristic function defining the region of integration
- class soupy.modeling.controlQoI.VariationalControlQoI(Vh, form_handler)[source]
Bases:
ControlQoIClass 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 variablesform_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 variablesrhs (
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 directiondir.- 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 variationout (
dolfin.Vector) – The assembled second variation
..note::
setLinearizationPointmust 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
ith,jth, andkth variables in directionsdir1anddir2.- Parameters:
i (int) – First variable index (
soupy.STATE,soupy.PARAMETER, orsoupy.CONTROL)j (int) – Second variable index (
soupy.STATE,soupy.PARAMETER, orsoupy.CONTROL)k (int) – Third variable index (
soupy.STATE,soupy.PARAMETER, orsoupy.CONTROL)dir1 (
dolfin.Vector) – Direction for variablejdir2 (
dolfin.Vector) – Direction for variablekout (
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
ith variable whereiis eithersoupy.STATE,soupy.PARAMETER, orsoupy.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 variablesout (
dolfin.Vector) – The assembled first variation
soupy.modeling.controlModelHessian
- class soupy.modeling.controlModelHessian.ControlModelHessian(model)[source]
Bases:
objectHessian of the the QoI map with respect to the
CONTROLvariable \(z\) wrapping around asoupy.ControlModelConstructor:
- Parameters:
model (
soupy.ControlModel) – Control model defining the control to QoI map
soupy.modeling.controlCostFunctional
- class soupy.modeling.controlCostFunctional.ControlCostFunctional[source]
Bases:
objectBase 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
componentisSTATE,PARAMETER,ADJOINT, orCONTROL, return a vector corresponding to that function space. Ifcomponentis"ALL", Generate the list of vectorsx = [u,m,p,z]
- class soupy.modeling.controlCostFunctional.DeterministicControlCostFunctional(model, prior, penalization=None)[source]
Bases:
ControlCostFunctionalThis 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 thesoupy.PDEVariationalControlProblemandsoupy.ControlQoIprior (
hippylib.Prior) – The prior distribution for the random parameterpenalization (
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 variableorder (int) – Order of the derivatives needed.
0for cost,1for gradient,2for Hessian
- cost(z, order=0, **kwargs)[source]
Computes the value of the cost functional at given control \(z\)
- Parameters:
z (
dolfin.Vector) – the control variableorder (int) – Order of the derivatives needed after evaluation
0for cost,1for gradient,2for Hessian
- Returns:
Value of the cost functional
- generate_vector(component='ALL')[source]
If
componentisSTATE,PARAMETER,ADJOINT, orCONTROL, return a vector corresponding to that function space. Ifcomponentis"ALL", Generate the list of vectorsx = [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.costhas been called withorder >= 1
- class soupy.modeling.controlCostFunctional.PenalizationControlCostFunctional(Vh, penalization, augment_control=False)[source]
Bases:
ControlCostFunctionalCost functional that is defined by a
soupy.Penalizationobject. 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
Trueif optimization problem uses asoupy.AugmentedVectoras 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
componentisSTATE,PARAMETER,ADJOINT,CONTROLreturn a vector corresponding to that function space. Ifcomponent == CONTROLandself.augment_controlisTrue, will return ansoupy.AugmentedVectorthat augments the control variablezwith a scalar that can be used for optimization. Ifcomponent == "ALL", Generate the list of vectorsx = [u,m,p,z]. Note that in this case, theCONTROLvariable will not be augmented with the scalar, and can be used directly for methods likesolveFwd.
- class soupy.modeling.controlCostFunctional.RiskMeasureControlCostFunctional(risk_measure, penalization=None)[source]
Bases:
objectThis 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 variableorder (int) – Order of the derivatives needed after evaluation
0for cost,1for gradient,2for Hessiankwargs – additional arguments, e.g.
rngfor the risk measure computation
- Returns:
Value of the cost functional
- generate_vector(component='ALL')[source]
If
componentisSTATE,PARAMETER,ADJOINT,CONTROLreturn a vector corresponding to that function space. Ifcomponent == CONTROLand the underlying risk measure uses asoupy.AugmentedVector, return ansoupy.AugmentedVectorthat augments the control variablezwith a scalar that can be used for optimization. Ifcomponent == "ALL", Generate the list of vectorsx = [u,m,p,z]. Note that in this case, theCONTROLvariable will not be augmented with the scalar, and can be used directly for methods likesolveFwd.
soupy.modeling.controlCostHessian
- class soupy.modeling.controlCostHessian.ControlCostHessian(cost_functional)[source]
Bases:
objectHessian of the the QoI map with respect to the
CONTROLvariable \(z\) wrapping around asoupy.ControlModel. Hessian will be applied wheresoupy.ControlModelis last called withorder>=2Constructor:
- 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.Vectoror similar) – Vector 1z2 (
dolfin.Vectoror 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.Vectoror similar) – Direction for Hessian actionHzhat (
dolfin.Vectoror similar) – Output for Hessian action
- property ncalls
soupy.modeling.riskMeasure
- class soupy.modeling.riskMeasure.RiskMeasure[source]
Bases:
objectAbstract 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 variableorder (int) – Order of the derivatives needed.
0for cost,1for gradient,2for Hessian
- cost()[source]
Evaluates the value of the risk measure once components have been computed
- Returns:
Value of the cost functional
Note
Assumes
computeComponentshas been called withorder>=0
- generate_vector(components='ALL')[source]
If
componentsisSTATE,PARAMETER,ADJOINT, orCONTROL, return a vector corresponding to that function space. Ifcomponentsis"ALL", Generate the list of vectorsx = [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.computeComponentshas been called withorder >= 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 measureHzhat (
dolfin.Vector) – (Dual of) Result of the Hessian action of the risk measure to store the result in
Note
Assumes
self.computeComponentshas been called withorder >= 1
soupy.modeling.meanVarRiskMeasureStochastic
- class soupy.modeling.meanVarRiskMeasureStochastic.MeanVarRiskMeasureStochastic(control_model, prior, settings=hippylib.ParameterList)[source]
Bases:
RiskMeasureClass 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 thesoupy.PDEVariationalControlProblemandsoupy.ControlQoIprior (
hippylib.Prior) – The prior distribution for the random parametersettings (
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 variableorder (int) – Order of the derivatives needed.
0for cost,1for gradient,2for Hessianrng (
hippylib.Random) – The random number generator used for sampling. IfNonethe default useshippylib.parRandom
- cost()[source]
Evaluates the value of the risk measure once components have been computed
- Returns:
Value of the cost functional
Note
Assumes
computeComponentshas been called withorder>=0
- generate_vector(component='ALL')[source]
If
componentsisSTATE,PARAMETER,ADJOINT, orCONTROL, return a vector corresponding to that function space. Ifcomponentsis"ALL", Generate the list of vectorsx = [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.computeComponentshas been called withorder >= 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 measureHzhat (
dolfin.Vector) – (Dual of) Result of the Hessian action of the risk measure to store the result in
Note
Assumes
self.computeComponentshas been called withorder >= 1
soupy.modeling.meanVarRiskMeasureSAA
- class soupy.modeling.meanVarRiskMeasureSAA.MeanVarRiskMeasureSAA(control_model, prior, settings=hippylib.ParameterList, comm_sampler=mpi4py.MPI.COMM_WORLD)[source]
Bases:
RiskMeasureMean + 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 thesoupy.PDEVariationalControlProblemandsoupy.ControlQoIprior (
hippylib.Prior) – The prior distribution for the random parametersettings (
hippylib.ParameterList) – additional settingscomm_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 variableorder (int) – Order of the derivatives needed.
0for cost,1for gradient,2for Hessiankwargs – 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
computeComponentshas been called withorder>=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
componentsisSTATE,PARAMETER,ADJOINT, orCONTROL, return a vector corresponding to that function space. Ifcomponentsis"ALL", generate the list of vectorsx = [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.computeComponentshas been called withorder >= 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 measureHzhat (
dolfin.Vector) – (Dual of) Result of the Hessian action of the risk measure to store the result in
Note
Assumes
self.computeComponentshas been called withorder >= 1
soupy.modeling.augmentedVector
- class soupy.modeling.augmentedVector.AugmentedVector(v, copy_vector=True)[source]
Bases:
objectClass 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 thetvariableNote
Currently supports only sample parallelism (i.e. serial mesh)
- Parameters:
v –
dolfin.Vectorto be augmentedcopy_vector – If
True, copy the vector, otherwise use the same memory
soupy.modeling.superquantileRiskMeasureSAA
- class soupy.modeling.superquantileRiskMeasureSAA.SuperquantileRiskMeasureSAA(control_model, prior, settings=hippylib.ParameterList, comm_sampler=mpi4py.MPI.COMM_WORLD)[source]
Bases:
RiskMeasureRisk 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 thesoupy.PDEVariationalControlProblemandsoupy.ControlQoIprior (
hippylib.Prior) – The prior distribution for the random parametersettings (
hippylib.ParameterList) – additional settings given insuperquantileRiskMeasureSAASettings()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 scalartappendedorder (int) – Order of the derivatives needed.
0for cost,1for gradient,2for Hessiankwargs – 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
computeComponentshas been called withorder>=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
componentisSTATE,PARAMETER,ADJOINT, return a vector corresponding to that function space. Ifcomponent == CONTROL, return ansoupy.AugmentedVectorthat augments the control variablezwith a scalar that can be used for optimizationIf
component == "ALL", Generate the list of vectorsx = [u,m,p,z]. Note that in this case, theCONTROLvariable will not be augmented with the scalar, and can be used directly for methods likesolveFwd.
- 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.computeComponentshas been called withorder>=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 measureHzt_hat (
soupy.AugmentedVector) – (Dual of) Result of the Hessian action of the risk measure to store the result in
Note
Assumes
self.computeComponentshas been called withorder >= 1
- 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:
objectWrapper 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
- class soupy.modeling.transformedMeanRiskMeasureSAA.IdentityFunction[source]
Bases:
objectScalar function \(f(x) = x\). Supports both
xas a single scalar or a vector of scalars.
- class soupy.modeling.transformedMeanRiskMeasureSAA.TransformedMeanRiskMeasureSAA(control_model, prior, settings=hippylib.ParameterList, comm_sampler=mpi4py.MPI.COMM_WORLD)[source]
Bases:
RiskMeasureRisk 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 thesoupy.PDEVariationalControlProblemandsoupy.ControlQoIprior (
hippylib.Prior) – The prior distribution for the random parametersettings (
hippylib.ParameterList) – Additional settingscomm_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 variableorder (int) – Order of the derivatives needed.
0for cost,1for gradient,2for Hessiankwargs – 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
computeComponentshas been called withorder>=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
componentisSTATE,PARAMETER,ADJOINT, orCONTROL, return a vector corresponding to that function space. Ifcomponentis"ALL", Generate the list of vectorsx = [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.computeComponentshas been called withorder >= 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 measureHzhat (
dolfin.Vector) – (Dual of) Result of the Hessian action of the risk measure to store the result in
Note
Assumes
self.computeComponentshas been called withorder >= 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 adolfin.VectorFunctionSpaceof realsConstructor
- Parameters:
Vh (list of
dolfin.FunctionSpace) – List of function spaces the state, parameter, adjoint, and controlalpha (float) – Weight for the penalization term
- cost(z)[source]
Compute the penalization at given control \(z\)
- Parameters:
z (
dolfin.VectororAugmentedVector) – The control variable
- grad(z, out)[source]
Compute the gradient of the penalty at given control \(z\)
- Parameters:
z (
dolfin.VectororAugmentedVector) – The control variableout (
dolfin.VectororAugmentedVector) – The assembled gradient vector
- hessian(z, zhat, out)[source]
Compute the Hessian of the penalty at given control \(z\)
- Parameters:
z (
dolfin.Vectororsoupy.AugmentedVector) – The control variablezhat (
dolfin.Vectororsoupy.AugmentedVector) – The direction for Hessian actionout (
dolfin.Vectororsoupy.AugmentedVector) – The assembled Hessian action vector
- class soupy.modeling.penalization.MultiPenalization(Vh, penalization_list, alpha_list=None)[source]
Bases:
PenalizationPenalization 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 controlpenalization_list (list of
Penalization) – List of penalization objectsalpha_list (list of floats) – List of weights for each penalization term
- cost(z)[source]
Compute the penalization at given control \(z\)
- Parameters:
z (
dolfin.Vectororsoupy.AugmentedVector) – The control variable
- grad(z, out)[source]
Compute the gradient of the penalty at given control \(z\)
- Parameters:
z (
dolfin.Vectororsoupy.AugmentedVector) – The control variableout (
dolfin.Vectororsoupy.AugmentedVector) – The assembled gradient vector
- hessian(z, zhat, out)[source]
Compute the Hessian of the penalty at given control \(z\)
- Parameters:
z (
dolfin.Vectororsoupy.AugmentedVector) – The control variablezhat (
dolfin.Vectororsoupy.AugmentedVector) – The direction for Hessian actionout (
dolfin.Vectororsoupy.AugmentedVector) – The assembled Hessian action vector
- class soupy.modeling.penalization.Penalization[source]
Bases:
objectAbstract class for the penalization on the control variable \(z\)
- class soupy.modeling.penalization.VariationalPenalization(Vh, form_handler)[source]
Bases:
PenalizationPenalization 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 variablesform_handler – The form handler for the penalization with a
__call__method that takes the control variablezas adolfin.Functionand returns the variational form.
- 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 controlM (
dolfin.Matrixlike) – Weighting matrix with methodmultalpha (float) – Weight for the penalization term
- cost(z)[source]
Compute the penalization at given control \(z\)
- Parameters:
z (
dolfin.Vectororsoupy.AugmentedVector) – The control variable
- grad(z, out)[source]
Compute the gradient of the penalty at given control \(z\)
- Parameters:
z (
dolfin.Vectororsoupy.AugmentedVector) – The control variableout (
dolfin.Vectororsoupy.AugmentedVector) – The assembled gradient vector
- hessian(z, zhat, out)[source]
Compute the Hessian of the penalty at given control \(z\)
- Parameters:
z (
dolfin.Vectororsoupy.AugmentedVector) – The control variablezhat (
dolfin.Vectororsoupy.AugmentedVector) – The direction for Hessian actionout (
dolfin.Vectororsoupy.AugmentedVector) – The assembled Hessian action vector
soupy.modeling.variables
- Enumerator for the variables of the inverse problem:
STATE = 0for the state variablePARAMETER = 1for the parameter variableADJOINT = 2for the adjoint variableCONTROL = 3for the control variable