ChainRules

PDEConstrainedFunctionals

GridapTopOpt.PDEConstrainedFunctionalsType
struct PDEConstrainedFunctionals{N,A}

An object that computes the objective, constraints, and their derivatives.

Implementation

This implementation computes derivatives of a integral quantity

$F(u(\varphi),\varphi,\mathrm{d}\Omega_1,\mathrm{d}\Omega_2,...) = \Sigma_{i}\int_{\Omega_i} f_i(\varphi)~\mathrm{d}\Omega$

with respect to an auxiliary parameter $\varphi$ where $u$ is the solution to a PDE and implicitly depends on $\varphi$. This requires two pieces of information:

  1. Computation of $\frac{\partial F}{\partial u}$ and $\frac{\partial F}{\partial \varphi}$ (handled by StateParamIntegrandWithMeasure).

  2. Computation of $\frac{\partial F}{\partial u} \frac{\partial u}{\partial \varphi}$ at $\varphi$ and $u$ using the adjoint method (handled by AbstractFEStateMap). I.e., let

    $\frac{\partial F}{\partial u} \frac{\partial u}{\partial \varphi} = -\lambda^\intercal \frac{\partial \mathcal{R}}{\partial \varphi}$

    where $\mathcal{R}$ is the residual and solve the (linear) adjoint problem:

    $\frac{\partial \mathcal{R}}{\partial u}^\intercal\lambda = \frac{\partial F}{\partial u}^\intercal.$

The gradient is then $\frac{\partial F}{\partial \varphi} = \frac{\partial F}{\partial \varphi} - \frac{\partial F}{\partial u}\frac{\partial u}{\partial \varphi}$.

Parameters

  • J: A StateParamIntegrandWithMeasure corresponding to the objective.
  • C: A vector of StateParamIntegrandWithMeasure corresponding to the constraints.
  • dJ: The DoFs for the objective sensitivity.
  • dC: The DoFs for each constraint sensitivity.
  • analytic_dJ: a Function for computing the analytic objective sensitivity.
  • analytic_dC: A vector of Function for computing the analytic objective sensitivities.
  • state_map::A: The state map for the problem.

Note

  • If analytic_dJ = nothing automatic differentiation will be used.
  • If analytic_dC[i] = nothing automatic differentiation will be used for C[i].
source
Gridap.Arrays.evaluate!Function
Fields.evaluate!(pcf::PDEConstrainedFunctionals,φh)

Evaluate the objective and constraints, and their derivatives at φh.

source

StateParamIntegrandWithMeasure

GridapTopOpt.StateParamIntegrandWithMeasureType
struct StateParamIntegrandWithMeasure{A<:IntegrandWithMeasure,B,C,D}

A wrapper to handle partial differentation of an IntegrandWithMeasure of a specific form (see below) in a ChainRules.jl compatible way with caching.

Assumptions

We assume that we have a IntegrandWithMeasure of the following form:

F(u,φ,dΩ₁,dΩ₂,...) = ∫(f(u,φ))dΩ₁ + ∫(g(u,φ))dΩ₂ + ...,.

where u and φ are each expected to inherit from Union{FEFunction,MultiFieldFEFunction} or the GridapDistributed equivalent.

source
ChainRulesCore.rruleMethod
ChainRulesCore.rrule(u_to_j::StateParamIntegrandWithMeasure,uh,φh)

Return the evaluation of a StateParamIntegrandWithMeasure and a a function for evaluating the pullback of u_to_j. This enables compatiblity with ChainRules.jl

source

Implemented types of AbstractFEStateMap

GridapTopOpt.AbstractFEStateMapType
abstract type AbstractFEStateMap

Types inheriting from this abstract type should enable the evaluation and differentiation of the solution to an FE problem u that implicitly depends on an auxiliary parameter φ.

source

AffineFEStateMap

GridapTopOpt.AffineFEStateMapType
struct AffineFEStateMap{A,B,C,D,E,F} <: AbstractFEStateMap

A structure to enable the forward problem and pullback for affine finite element operators AffineFEOperator.

Parameters

  • biform::A: IntegrandWithMeasure defining the bilinear form.
  • liform::B: IntegrandWithMeasure defining the linear form.
  • spaces::C: Tuple of finite element spaces.
  • plb_caches::D: A cache for the pullback operator.
  • fwd_caches::E: A cache for the forward problem.
  • adj_caches::F: A cache for the adjoint problem.
source
GridapTopOpt.AffineFEStateMapMethod
AffineFEStateMap(
  a::Function,l::Function,
  U,V,V_φ,U_reg,φh,dΩ...;
  assem_U = SparseMatrixAssembler(U,V),
  assem_adjoint = SparseMatrixAssembler(V,U),
  assem_deriv = SparseMatrixAssembler(U_reg,U_reg),
  ls::LinearSolver = LUSolver(),
  adjoint_ls::LinearSolver = LUSolver()
)

Create an instance of AffineFEStateMap given the bilinear form a and linear form l as Function types, trial and test spaces U and V, the FE space V_φ for φh, the FE space U_reg for derivatives, and the measures as additional arguments.

Optional arguments enable specification of assemblers and linear solvers.

source

NonlinearFEStateMap

GridapTopOpt.NonlinearFEStateMapType
struct NonlinearFEStateMap{A,B,C,D,E,F} <: AbstractFEStateMap

A structure to enable the forward problem and pullback for nonlinear finite element operators.

Parameters

  • res::A: an IntegrandWithMeasure defining the residual of the problem.
  • jac::B: a Function defining Jacobian of the residual.
  • spaces::C: Tuple of finite element spaces.
  • plb_caches::D: A cache for the pullback operator.
  • fwd_caches::E: A cache for the forward problem.
  • adj_caches::F: A cache for the adjoint problem.
source
GridapTopOpt.NonlinearFEStateMapMethod
NonlinearFEStateMap(
  res::Function,U,V,V_φ,U_reg,φh,dΩ...;
  assem_U = SparseMatrixAssembler(U,V),
  assem_adjoint = SparseMatrixAssembler(V,U),
  assem_deriv = SparseMatrixAssembler(U_reg,U_reg),
  nls::NonlinearSolver = NewtonSolver(LUSolver();maxiter=50,rtol=1.e-8,verbose=true),
  adjoint_ls::LinearSolver = LUSolver()
)

Create an instance of NonlinearFEStateMap given the residual res as a Function type, trial and test spaces U and V, the FE space V_φ for φh, the FE space U_reg for derivatives, and the measures as additional arguments.

Optional arguments enable specification of assemblers, nonlinear solver, and adjoint (linear) solver.

source

RepeatingAffineFEStateMap

GridapTopOpt.RepeatingAffineFEStateMapType
struct RepeatingAffineFEStateMap <: AbstractFEStateMap

A structure to enable the forward problem and pullback for affine finite element operators AffineFEOperator with multiple linear forms but only a single bilinear form.

Parameters

  • biform: IntegrandWithMeasure defining the bilinear form.
  • liform: A vector of IntegrandWithMeasure defining the linear forms.
  • spaces: Repeated finite element spaces.
  • spaces_0: Original finite element spaces that are being repeated.
  • plb_caches: A cache for the pullback operator.
  • fwd_caches: A cache for the forward problem.
  • adj_caches: A cache for the adjoint problem.
source
GridapTopOpt.RepeatingAffineFEStateMapMethod
RepeatingAffineFEStateMap(
  nblocks::Int,a::Function,l::Vector{<:Function},
  U0,V0,V_φ,U_reg,φh,dΩ...;
  assem_U = SparseMatrixAssembler(U0,V0),
  assem_adjoint = SparseMatrixAssembler(V0,U0),
  assem_deriv = SparseMatrixAssembler(U_reg,U_reg),
  ls::LinearSolver = LUSolver(),
  adjoint_ls::LinearSolver = LUSolver()
)

Create an instance of RepeatingAffineFEStateMap given the number of blocks nblocks, a bilinear form a, a vector of linear form l as Function types, the trial and test spaces U and V, the FE space V_φ for φh, the FE space U_reg for derivatives, and the measures as additional arguments.

Optional arguments enable specification of assemblers and linear solvers.

Note

  • The resulting FEFunction will be a MultiFieldFEFunction (or GridapDistributed equivalent) where each field corresponds to an entry in the vector of linear forms
source

Advanced

Inheriting from AbstractFEStateMap

Existing methods

ChainRulesCore.rruleMethod
rrule(φ_to_u::AbstractFEStateMap,φh)

Return the evaluation of a AbstractFEStateMap and a a function for evaluating the pullback of φ_to_u. This enables compatiblity with ChainRules.jl

source
GridapTopOpt.pullbackFunction
pullback(φ_to_u::AbstractFEStateMap,uh,φh,du;updated)

Compute ∂F∂u*dudφ at φh and uh using the adjoint method. I.e., let

∂F∂u*dudφ = -λᵀ*dRdφ

and solve the adjoint problem dRduᵀ*λ = ∂F∂uᵀ using adjoint_solve!.

source

Required to implement

GridapTopOpt.forward_solve!Function
forward_solve!(φ_to_u::AbstractFEStateMap,φh)

Evaluate the forward problem u given φ. This should compute the FE problem.

source
GridapTopOpt.adjoint_solve!Function
adjoint_solve!(φ_to_u::AbstractFEStateMap,du::AbstractVector)

Evaluate the solution to the adjoint problem given a RHS vector ∂F∂u denoted du. This should solve the linear problem dRduᵀ*λ = ∂F∂uᵀ.

source
GridapTopOpt.dRdφFunction
dRdφ(φ_to_u::AbstractFEStateMap,uh,vh,φh)

Compute the derivative with respect to φh of the residual R.

source

IntegrandWithMeasure

GridapTopOpt.IntegrandWithMeasureType
struct IntegrandWithMeasure{A,B<:Tuple}

A wrapper to enable serial or parallel partial differentation of an integral F using Gridap.gradient. This is required to allow automatic differentation with DistributedMeasure.

Properties

  • F :: A: A function that returns a DomainContribution or DistributedDomainContribution.
  • dΩ :: B: A tuple of measures.
source
Gridap.Fields.gradientFunction
Gridap.gradient(F::IntegrandWithMeasure,uh::Vector,K::Int)

Given an an IntegrandWithMeasure F and a vector of FEFunctions uh (excluding measures) evaluate the partial derivative of F.F with respect to uh[K].

Example

Suppose uh and φh are FEFunctions with measures and dΓ_N. Then the partial derivative of a function J wrt to φh is computed via

J(u,φ,dΩ,dΓ_N) = ∫(f(u,φ))dΩ + ∫(g(u,φ))dΓ_N
J_iwm = IntegrandWithMeasure(J,(dΩ,dΓ_N))
∂J∂φh = ∇(J_iwm,[uh,φh],2)

where f and g are user defined.

source
Gridap.Algebra.jacobianFunction
Gridap.jacobian(F::IntegrandWithMeasure,uh::Vector,K::Int)

Given an an IntegrandWithMeasure F and a vector of FEFunctions or CellField uh (excluding measures) evaluate the Jacobian F.F with respect to uh[K].

source