gemseo / core

jacobian_assembly module

Coupled derivatives calculations.

Classes:

CoupledSystem()

Compute coupled (total) derivatives of a system of residuals.

JacobianAssembly(coupling_structure)

Assembly of Jacobians.

Functions:

default_dict_factory()

Instantiates a defaultdict(None) object.

none_factory()

Returns None...

class gemseo.core.jacobian_assembly.CoupledSystem[source]

Bases: object

Compute coupled (total) derivatives of a system of residuals.

Use several methods:

  • direct or adjoint

  • factorized for multiple RHS

Methods:

adjoint_mode(functions, dres_dx, dres_dy_t, ...)

Compute the total derivative Jacobian in adjoint mode.

direct_mode(functions, n_variables, ...[, ...])

Computate the total derivative Jacobian in direct mode.

adjoint_mode(functions, dres_dx, dres_dy_t, dfun_dx, dfun_dy, linear_solver='DEFAULT', use_lu_fact=False, **linear_solver_options)[source]

Compute the total derivative Jacobian in adjoint mode.

Parameters
  • functions – The functions to differentiate.

  • dres_dx – The Jacobian of the residuals wrt the design variables.

  • dres_dy – The Jacobian of the residuals wrt the coupling variables.

  • dfun_dx – The Jacobian of the functions wrt the design variables.

  • dfun_dy – The Jacobian of the functions wrt the coupling variables.

  • linear_solver

    The name of the linear solver.

    By default it is set to DEFAULT.

  • use_lu_fact

    Whether to factorize dres_dy_t once.

    By default it is set to False.

  • dres_dy_t – The param kwargs.

  • linear_solver_options – The optional parameters.

Returns

The Jacobian of total coupled derivatives.

direct_mode(functions, n_variables, n_couplings, dres_dx, dres_dy, dfun_dx, dfun_dy, linear_solver='DEFAULT', use_lu_fact=False, **linear_solver_options)[source]

Computate the total derivative Jacobian in direct mode.

Parameters
  • functions – The functions to differentiate.

  • n_variables – The number of variables.

  • n_couplings – The number of couplings.

  • dres_dx – The Jacobian of the residuals wrt the design variables.

  • dres_dy – The Jacobian of the residuals wrt the coupling variables.

  • dfun_dx – The Jacobian of the functions wrt the design variables.

  • dfun_dy – The Jacobian of the functions wrt the coupling variables.

  • linear_solver

    The name of the linear solver.

    By default it is set to DEFAULT.

  • use_lu_fact

    Whether to factorize dres_dy once.

    By default it is set to False.

  • linear_solver_options – The optional parameters.

Returns

The Jacobian of the total coupled derivatives.

class gemseo.core.jacobian_assembly.JacobianAssembly(coupling_structure)[source]

Bases: object

Assembly of Jacobians.

Typically, assemble disciplines’s Jacobians into a system Jacobian.

Parameters
  • coupling_structure – The CouplingStructure associated disciplines that form

  • system. (the coupled) –

Attributes:

ADJOINT_MODE

AUTO_MODE

AVAILABLE_MAT_TYPES

AVAILABLE_MODES

DIRECT_MODE

LINEAR_OPERATOR

REVERSE_MODE

SPARSE

Methods:

compute_dimension(names)

Compute the total number of functions/variables/couplings of the whole system.

compute_newton_step(in_data, couplings, ...)

Compute the Newton step for the coupled system of residuals formed by the disciplines.

compute_sizes(functions, variables, couplings)

Compute the number of scalar functions, variables and couplings.

dfun_dvar(function, variables, n_variables)

Forms the matrix of partial derivatives of a function.

dres_dvar(residuals, variables, n_residuals, ...)

Form the matrix of partial derivatives of residuals.

plot_dependency_jacobian(functions, variables)

Plot the Jacobian matrix Nonzero elements of the sparse matrix are represented by blue squares.

residuals(in_data, var_names)

Form the matrix of residuals wrt coupling variables.

split_jac(coupled_system, variables)

Split a Jacobian dict into a dict of dict.

total_derivatives(in_data, functions, ...[, ...])

Compute the Jacobian of total derivatives of the coupled system formed by the disciplines.

ADJOINT_MODE = 'adjoint'
AUTO_MODE = 'auto'
AVAILABLE_MAT_TYPES = ['sparse', 'linear_operator']
AVAILABLE_MODES = ('direct', 'adjoint', 'auto', 'reverse')
DIRECT_MODE = 'direct'
LINEAR_OPERATOR = 'linear_operator'
REVERSE_MODE = 'reverse'
SPARSE = 'sparse'
compute_dimension(names)[source]

Compute the total number of functions/variables/couplings of the whole system.

Parameters

names – The names of the inputs or the outputs.

Returns

The dimension if the system.

compute_newton_step(in_data, couplings, relax_factor, linear_solver='DEFAULT', matrix_type='sparse', **linear_solver_options)[source]

Compute the Newton step for the coupled system of residuals formed by the disciplines.

Parameters
  • in_data – The input data.

  • couplings – The coupling variables.

  • relax_factor – The relaxation factor.

  • linear_solver

    The name of the linear solver.

    By default it is set to DEFAULT.

  • matrix_type

    The representation of the matrix dR/dy (sparse or linear operator).

    By default it is set to sparse.

  • linear_solver_options – The options passed to the linear solver factory.

Returns

The Newton step -[dR/dy]^-1 . R as a dict of steps per coupling variable.

compute_sizes(functions, variables, couplings)[source]

Compute the number of scalar functions, variables and couplings.

Parameters
  • functions – The functions to differentiate.

  • variables – The differentiation variables.

  • couplings – The coupling variables.

Raises

ValueError – When the size of some variables could not be determined.

dfun_dvar(function, variables, n_variables)[source]

Forms the matrix of partial derivatives of a function.

Given disciplinary Jacobians dJi(v0…vn)/dvj, fill the sparse Jacobian: | | | dJi/dvj | | |

Parameters
  • function – The function to differentiate.

  • variables – The differentiation variables.

  • n_variables – The number of variables.

Returns

The full Jacobian matrix.

dres_dvar(residuals, variables, n_residuals, n_variables, matrix_type='sparse', transpose=False)[source]

Form the matrix of partial derivatives of residuals.

Given disciplinary Jacobians dYi(Y0…Yn)/dvj, fill the sparse Jacobian: | | | dRi/dvj | | | (Default value = False)

Parameters
  • residuals – The residuals.

  • variables – The differentiation variables.

  • n_residuals – The number of residuals.

  • n_variables – The number of variables.

  • matrix_type

    The type of the matrix.

    By default it is set to sparse.

  • transpose

    Whether to transpose the matrix.

    By default it is set to False.

Returns

The jacobian of dres_dvar.

Raises

TypeError – When the matrix type is unknown.

plot_dependency_jacobian(functions, variables, save=True, show=False, filepath=None, markersize=None)[source]

Plot the Jacobian matrix Nonzero elements of the sparse matrix are represented by blue squares.

Parameters
  • functions – The functions to plot.

  • variables – The variables.

  • show

    WHether the plot is displayed.

    By default it is set to False.

  • save

    WHether the plot is saved in a PDF file.

    By default it is set to True.

  • filepath

    The file name to save to. If None, “coupled_jacobian.pdf” is used, otherwise “coupled_jacobian_” + filepath + “.pdf”.

    By default it is set to None.

Returns

The file name.

residuals(in_data, var_names)[source]

Form the matrix of residuals wrt coupling variables.

Given disciplinary explicit calculations Yi(Y0_t,…Yn_t), fill the residual matrix:

[Y0(Y0_t,...Yn_t) - Y0_t]
[                       ]
[Yn(Y0_t,...Yn_t) - Yn_t]
Parameters
  • in_data – The values prescribed for the calculation of the residuals (Y0_t,…Yn_t).

  • var_names – The names of variables associated with the residuals (R).

Returns

The residuals array.

split_jac(coupled_system, variables)[source]

Split a Jacobian dict into a dict of dict.

Parameters
  • coupled_system – The derivatives to split.

  • variables – The variables wrt which the differentiation is performed.

Returns

The Jacobian.

total_derivatives(in_data, functions, variables, couplings, linear_solver='DEFAULT', mode='auto', matrix_type='sparse', use_lu_fact=False, exec_cache_tol=None, force_no_exec=False, **linear_solver_options)[source]

Compute the Jacobian of total derivatives of the coupled system formed by the disciplines.

Parameters
  • in_data – The input data dict.

  • functions – The functions to differentiate.

  • variables – The differentiation variables.

  • couplings – The coupling variables.

  • linear_solver

    The name of the linear solver.

    By default it is set to DEFAULT.

  • mode

    The linearization mode (auto, direct or adjoint).

    By default it is set to auto.

  • matrix_type

    The representation of the matrix dR/dy (sparse or linear operator).

    By default it is set to sparse.

  • use_lu_fact

    Whether to factorize dres_dy once, unsupported for linear operator mode.

    By default it is set to False.

  • exec_cache_tol

    The discipline cache tolerance to when calling the linearize method. If None, no tolerance is set (equivalent to tol=0.0).

    By default it is set to None.

  • force_no_exec

    Whether the discipline is not re executed, the cache is loaded anyway.

    By default it is set to False.

  • linear_solver_options – The options passed to the linear solver factory.

Returns

The total coupled derivatives.

Raises

ValueError – When the linearization_mode is incorrect.

gemseo.core.jacobian_assembly.default_dict_factory()[source]

Instantiates a defaultdict(None) object.

gemseo.core.jacobian_assembly.none_factory()[source]

Returns None…

To be used for defaultdict