Skip to content

Adjoint API

Michael Andre edited this page Jul 12, 2017 · 24 revisions

Overview

This page describes the API for adjoint-based sensitivities in Kratos. It is the reference for developers creating new adjoint elements and objective functions.

Adjoint Variables

An adjoint variable is defined for each primal variable by prefixing the variable name with ADJOINT. For example, the adjoint variable of DISPLACEMENT is ADJOINT_DISPLACEMENT.

Adjoint Elements and Conditions

Adjoint elements are created in the same way as regular elements. The adjoint element should be a separate class from the primal element. For applications with a relatively large number of linear elements like structural mechanics, the primal element type can be given as a template argument to a single generic adjoint element which stores the element internally and uses the element's public member functions to produce the adjoint matrices. Inheritance can also be used to minimize code duplication between primal and adjoint elements. The element functions and their meaning for the adjoint problem are given below.

Solving the adjoint problem

The following functions are used in the solution of the adjoint problem:

GetDofList(DofsVectorType& rElementalDofList, ProcessInfo& rCurrentProcessInfo)

Gets the list of adjoint dofs. This is the analog of the list of primal dofs.

CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix, ProcessInfo& rCurrentProcessInfo)

Calculates adjoint matrix as the gradient of the negative residual with respect to the primal variable (e.g., DISPLACEMENT). For structural mechanics, this is the transpose of the stiffness matrix. This is zero if the primal variable is a rate (e.g., VELOCITY for fluid mechanics).

CalculateFirstDerivativesLHS(MatrixType& rLeftHandSideMatrix, ProcessInfo& rCurrentProcessInfo)

Calculates the adjoint matrix as the gradient of the negative residual with respect to first derivatives of the primal variable (e.g., VELOCITY). For structural mechanics, this is the transpose of the damping matrix.

CalculateSecondDerivativesLHS(MatrixType& rLeftHandSideMatrix, ProcessInfo& rCurrentProcessInfo)

Calculates the adjoint matrix as the gradient of the negative residual with respect to second derivatives of the primal variable (e.g., ACCELERATION). This is normally the transpose of the mass matrix.

GetValuesVector(Vector& values, int Step = 0)

Gets the vector of adjoint values for the primal variable (e.g., ADJOINT_DISPLACEMENT).

GetFirstDerivativesVector(Vector& values, int Step = 0)

Gets the vector of adjoint values for first derivatives of primal variable (e.g., ADJOINT_VELOCITY).

GetSecondDerivativesVector(Vector& values, int Step = 0)

Gets the vector of adjoint values for second derivatives of primal variable (e.g., ADJOINT_ACCELERATION).

Calculating sensitivities

The following functions are used after the solution of the adjoint problem to calculate sensitivities:

CalculateSensitivityMatrix(const Variable<double>& rVariable, Matrix& Output, const ProcessInfo& rCurrentProcessInfo)
CalculateSensitivityMatrix(const Variable<array_1d<double,3>>& rVariable, Matrix& Output, const ProcessInfo& rCurrentProcessInfo)

Calculates the transpose of the gradient of the element residual with respect to the design variable. (These functions must be added to the element base class).

The sensitivities depend on the design variable which may be a double (e.g., THICKNESS) or array_1d<double,3> (e.g., NODAL_COORDINATES). The row dimension of the output is determined from the dimension of design variable (scalar or array_1d) and its storage type (nodal or elemental).

Response Functions

The response function defines a scalar quantity of interest (e.g., drag), its gradients with respect to primal variables and their derivatives and its gradient with respect to design variables.

The following functions should be defined for each response function:

CalculateValue()

Calculates the scalar value of the response function.

CalculateGradient(const Element& rElem, const Matrix& rLHS, Vector& rOutput, ProcessInfo& rProcessInfo)

Calculates the gradient of the response function with respect to the primal variable (e.g., DISPLACEMENT).

CalculateFirstDerivativesGradient(const Element& rElem, const Matrix& rLHS, Vector& rOutput, ProcessInfo& rProcessInfo)

Calculates the gradient of the response function with respect to first derivatives of the primal variable (e.g., VELOCITY).

CalculateSecondDerivativesGradient(const Element& rElem, const Matrix& rLHS, Vector& rOutput, ProcessInfo& rProcessInfo)

Calculates the gradient of the response function with respect to second derivatives of the primal variable (e.g., ACCELERATION).

CalculateSensitivityGradient(const Variable<double>& rVariable, const Element& rElem, const Matrix& rSensitivityMatrix, Vector& rOutput, ProcessInfo& rProcessInfo)
CalculateSensitivityGradient(const Variable<array_1d<double,3>>& rVariable, const Element& rElem, const Matrix& rSensitivityMatrix, Vector& rOutput, ProcessInfo& rProcessInfo)

Calculates the gradient of the response function with respect to design parameters.

Project information

Getting Started

Tutorials

Developers

Kratos structure

Conventions

Solvers

Debugging, profiling and testing

HOW TOs

Utilities

Kratos API

Kratos Structural Mechanics API

Clone this wiki locally