Reference

Contents

Index

PDENLPModels.EnergyFETermType

AbstractEnergyTerm modeling the objective function of the optimization problem

\[\begin{aligned} \int_{\Omega} f(y,u) d\Omega. \end{aligned}\]

Constructor:

EnergyFETerm(:: Function, :: Triangulation, :: Measure)

See also: MixedEnergyFETerm, NoFETerm, _obj_cell_integral, _obj_integral, _compute_gradient_k!

source
PDENLPModels.GridapPDENLPModelType

GridapPDENLPModel returns an instance of an AbstractNLPModel using Gridap.jl for the discretization of the domain with finite-elements. Given a domain Ω ⊂ ℜᵈ Find a state function y: Ω -> Y, a control function u: Ω -> U and an algebraic vector κ ∈ ℜⁿ satisfying

\[\begin{aligned} \min_{κ,y,u} \ & ∫_Ω​ f(κ,y,u) dΩ​ \\ \mbox{ s.t. } & y \mbox{ solution of } PDE(κ,u)=0, \\ & lcon <= c(κ,y,u) <= ucon, \\ & lvar <= (κ,y,u) <= uvar. \end{aligned}\]

The weak formulation of the PDE is then: res((y,u),(v,q)) = ∫ v PDE(κ,y,u) + ∫ q c(κ,y,u)

where the unknown (y,u) is a MultiField see Tutorials 7 and 8 of Gridap.

Constructors

Main constructor:

GridapPDENLPModel(::NLPModelMeta, ::Counters, ::PDENLPMeta, ::PDEWorkspace)

This is the main constructors with the attributes of the GridapPDENLPModel:

  • meta::NLPModelMeta: usual meta for NLPModels, see doc here;
  • counters::Counters: usual counters for NLPModels, see doc here;
  • pdemeta::PDENLPMeta: metadata specific to GridapPDENLPModel;
  • workspace::PDEWorkspace: Pre-allocated memory for GridapPDENLPModel.

More practical constructors are also available.

  • For unconstrained problems:

    GridapPDENLPModel(x0, ::Function, ::Union{Measure, Triangulation}, Ypde::FESpace, Xpde::FESpace; kwargs...)

    GridapPDENLPModel(x0, ::AbstractEnergyTerm, ::Union{Measure, Triangulation}, Ypde::FESpace, Xpde::FESpace; kwargs...)

  • For constrained problems without controls:

    GridapPDENLPModel(x0, ::Function, ::Union{Measure, Triangulation}, Ypde::FESpace, Xpde::FESpace, c::Union{Function, FEOperator}; kwargs...)

    GridapPDENLPModel(x0, ::AbstractEnergyTerm, Ypde::FESpace, Xpde::FESpace, c::Union{Function, FEOperator}; kwargs...)

  • For general constrained problems:

    GridapPDENLPModel(x0, ::Function, ::Union{Measure, Triangulation}, Ypde::FESpace, Ycon::FESpace, Xpde::FESpace, Xcon::FESpace, c::Union{Function, FEOperator}; kwargs...)

    GridapPDENLPModel(x0, ::AbstractEnergyTerm, Ypde::FESpace, Ycon::FESpace, Xpde::FESpace, Xcon::FESpace, c::Union{Function, FEOperator}; kwargs...)

where the different arguments are:

  • x0: initial guess for the system of size ≥ num_free_dofs(Ypde) + num_free_dofs(Ycon);
  • f: objective function, the number of arguments depend on the application (y) or (y,u) or (y,u,θ);
  • Ypde: trial space for the state;
  • Ycon: trial space for the control (VoidFESpace if none);
  • Xpde: test space for the state;
  • Xcon: test space for the control (VoidFESpace if none);
  • c: operator/function for the PDE-constraint, were we assume by default that the right-hand side is zero (otw. use lcon and ucon keywords), the number of arguments depend on the application (y,v) or (y,u,v) or (y,u,θ,v).

If length(x0) > num_free_dofs(Ypde) + num_free_dofs(Ycon), then the additional components are considered algebraic variables.

The function f and c must return integrals complying with Gridap's functions with a Measure/Triangulation given in the arguments of GridapPDENLPModel. Internally, the objective function f and the Measure/Triangulation are combined to instantiate an AbstractEnergyTerm.

The following keyword arguments are available to all constructors:

  • name: The name of the model (default: "Generic")

The following keyword arguments are available to the constructors for constrained problems:

  • lin: An array of indexes of the linear constraints (default: Int[])

The following keyword arguments are available to the constructors for constrained problems explictly giving lcon and ucon:

  • y0: An inital estimate to the Lagrangian multipliers (default: zeros)

The bounds on the variables are given as AbstractVector via keywords arguments as well:

  • either with lvar and uvar, or,
  • lvary, lvaru, lvark, and uvary, uvaru, uvark.

Notes:

  • We handle two types of FEOperator: AffineFEOperator, and FEOperatorFromWeakForm.
  • If lcon and ucon are not given, they are assumed zeros.
  • If the type can't be deduced from the argument, it is Float64.

Example

using Gridap, PDENLPModels

  # Definition of the domain
  n = 100
  domain = (-1, 1, -1, 1)
  partition = (n, n)
  model = CartesianDiscreteModel(domain, partition)

  # Definition of the spaces:
  valuetype = Float64
  reffe = ReferenceFE(lagrangian, valuetype, 2)
  Xpde = TestFESpace(model, reffe; conformity = :H1, dirichlet_tags = "boundary")
  y0(x) = 0.0
  Ypde = TrialFESpace(Xpde, y0)

  reffe_con = ReferenceFE(lagrangian, valuetype, 1)
  Xcon = TestFESpace(model, reffe_con; conformity = :H1)
  Ycon = TrialFESpace(Xcon)

  # Integration machinery
  trian = Triangulation(model)
  degree = 1
  dΩ = Measure(trian, degree)

  # Objective function:
  yd(x) = -x[1]^2
  α = 1e-2
  function f(y, u)
    ∫(0.5 * (yd - y) * (yd - y) + 0.5 * α * u * u) * dΩ
  end

  # Definition of the constraint operator
  ω = π - 1 / 8
  h(x) = -sin(ω * x[1]) * sin(ω * x[2])
  function res(y, u, v)
    ∫(∇(v) ⊙ ∇(y) - v * u - v * h) * dΩ
  end

  # initial guess
  npde = num_free_dofs(Ypde)
  ncon = num_free_dofs(Ycon)
  xin = zeros(npde + ncon)

  nlp = GridapPDENLPModel(xin, f, trian, Ypde, Ycon, Xpde, Xcon, res, name = "Control elastic membrane")

You can also check the tutorial Solve a PDE-constrained optimization problem on JSO's website, juliasmoothoptimizers.github.io.

We refer to the folder test/problems for more examples of problems of different types:

  • calculus of variations,
  • optimal control problem,
  • PDE-constrained problems,
  • mixed PDE-contrained problems with both function and algebraic unknowns.

An alternative is to visit the repository PDEOptimizationProblems that contains a collection of test problems.

Without objective function, the problem reduces to a classical PDE and we refer to Gridap tutorials for examples.

source
PDENLPModels.MixedEnergyFETermType

AbstractEnergyTerm modeling the objective function of the optimization problem with functional and discrete unknowns

\[\begin{aligned} \int_{\Omega} f(y,u,\kappa) d\Omega. \end{aligned}\]

Constructor:

MixedEnergyFETerm(:: Function, :: Triangulation, :: Int)

See also: EnergyFETerm, NoFETerm, _obj_cell_integral, _obj_integral, _compute_gradient_k!

source
PDENLPModels.NoFETermType

AbstractEnergyTerm modeling the objective function when there are no integral objective

\[\begin{aligned} f(\kappa). \end{aligned}\]

Constructors: NoFETerm() NoFETerm(:: Function)

See also: MixedEnergyFETerm, EnergyFETerm, _obj_cell_integral, _obj_integral, _compute_gradient_k!

source
PDENLPModels.PDENLPMetaType
PDENLPMeta

A composite type that represents the main features of the PDE-constrained optimization problem


PDENLPMeta contains the following attributes:

  • tnrj : structure representing the objective term
  • Ypde : TrialFESpace for the solution of the PDE
  • Ycon : TrialFESpace for the parameter
  • Xpde : TestFESpace for the solution of the PDE
  • Xcon : TestFESpace for the parameter
  • Y : concatenated TrialFESpace
  • X : concatenated TestFESpace
  • op : operator representing the PDE-constraint (nothing if no constraints)
  • nvar_pde :number of dofs in the solution functions
  • nvar_con : number of dofs in the control functions
  • nparam : number of real unknowns
  • nnzh_obj : number of nonzeros elements in the objective hessian
  • Hrows : store the structure for the hessian of the lagrangian
  • Hcols : store the structure for the hessian of the lagrangian
  • Jrows : store the structure for the hessian of the jacobian
  • Jcols : store the structure for the hessian of the jacobian
source
PDENLPModels._compute_gradient!Function

Return the gradient of the objective function and set it in place.

_compute_gradient!(:: AbstractVector, :: EnergyFETerm, :: AbstractVector, :: FEFunctionType, :: FESpace, :: FESpace)

See also: MixedEnergyFETerm, EnergyFETerm, NoFETerm, _obj_integral, _obj_cell_integral, _compute_hess_k_coo

source
PDENLPModels._compute_gradient_k!Function

Return the derivative of the objective function w.r.t. κ.

_compute_gradient_k!(::AbstractVector, :: AbstractEnergyTerm, :: FEFunctionType, :: AbstractVector)

See also: MixedEnergyFETerm, EnergyFETerm, NoFETerm, _obj_integral, _obj_cell_integral, _compute_hess_k_coo

source
PDENLPModels._compute_hess_k_vals!Function

Return the values of the hessian w.r.t. κ of the objective function.

_compute_hess_k_vals!(:: AbstractVector, :: AbstractNLPModel, :: AbstractEnergyTerm, :: AbstractVector, :: AbstractVector)

See also: MixedEnergyFETerm, EnergyFETerm, NoFETerm, _obj_integral, _obj_cell_integral, _compute_gradient_k

source
PDENLPModels._compute_hess_structureMethod
_compute_hess_structure(::AbstractEnergyTerm, Y, X, x0, nparam)
_compute_hess_structure::AbstractEnergyTerm, op, Y, Ypde, Ycon, X, Xpde, x0, nparam)

Return a triplet with the structure (rows, cols) and the number of non-zeros elements in the hessian w.r.t. y and u of the objective function.

The rows and cols returned by _compute_hess_structure_obj are already shifter by nparam.

source
PDENLPModels._functions_to_vectors!Method
_functions_to_vectors!(nini :: Int, nfields :: Int, trian :: Triangulation, lfunc :: Function, ufunc :: Function, cell_xm :: Gridap.Arrays.AppliedArray, Y :: FESpace, lvar :: AbstractVector, uvar :: AbstractVector)

Iterate for k = 1 to nfields and switch lfunc[k] and ufunc[k] to vectors, allocated in lvar and uvar in place starting from nini + 1. It returns nini + the number of allocations.

source
PDENLPModels._obj_integralFunction

Return the integral of the objective function

_obj_integral(:: AbstractEnergyTerm, :: FEFunctionType, :: AbstractVector)

See also: MixedEnergyFETerm, EnergyFETerm, NoFETerm, _obj_cell_integral, _compute_gradient_k, _compute_hess_k_coo

source
PDENLPModels._split_FEFunctionMethod
_split_FEFunction(:: AbstractVector,  :: FESpace, :: Union{FESpace, Nothing})

Split the vector x into two FEFunction corresponding to the solution y and the control u. Returns nothing for the control u if Ycon == nothing.

Do not verify the compatible length.

source
PDENLPModels._split_vectorMethod
_split_vector(:: AbstractVector,  :: FESpace, :: Union{FESpace, Nothing})

Split the vector x into three vectors: y, u, k. Returns nothing for the control u if Ycon == nothing.

Do not verify the compatible length.

source
PDENLPModels.bounds_functions_to_vectorsMethod
bounds_functions_to_vectors(Y :: MultiFieldFESpace, Ycon :: Union{FESpace, Nothing},  Ypde :: FESpace, trian :: Triangulation, lyfunc :: Union{Function, AbstractVector}, uyfunc :: Union{Function, AbstractVector}, lufunc :: Union{Function, AbstractVector}, uufunc :: Union{Function, AbstractVector})

Return the bounds lvar and uvar.

source