source_localization.adjoint_convection_diffusion_2D module

class source_localization.adjoint_convection_diffusion_2D.AdjointDiffusionConvectionReaction2DEquation(U, K, coeff_depot, msh, time_discretization='semi-implicit', tol_inversion=1e-14)[source]

Bases: object

Class containing the adjoint model of the 2D diffusion-convection-reaction PDE model :

\[\partial_tc^* + \nabla\cdot(\mathbf{K}^T\nabla c^*)+\nabla(\vec{u}c^*)-(\nabla.\vec{u})c^*-\tau_{loss}c^* = \left(\frac{dm}{dc}(c(s))\right)^*\cdot2\mathbf{R}^{-1}\left(m(c(s))-m_{obs}\right) ~\forall (x,y)\in\Omega~\forall t\in[0;T[\]

with the final and boundary conditions:

  • a null final condition \(c(x,y,t=T)=0~\forall (x,y)\in\Omega\),

  • a null diffusive flux \(\mathbf{K}^T\nabla c\cdot\vec{n}=0~\forall (x,y)\in\partial\Omega\),

  • a null outgoing convective flux \(\vec{u}c\cdot\vec{n}=0~\forall (x,y)\in\partial\Omega\cap \{(x,y)|\vec{u}(x,y,t)\cdot\vec{n}>0\}~\forall t\in]0;T]\) with \(\vec{n}\) the outgoing normal vector,

and its solvers.

msh

The geometry of the domain.

Type:

MeshRect2D

A_adjoint

The adjoint of the convection linear operator \(A^*:c^*\mapsto - \nabla\cdot(\vec{u}c^*)+(\nabla.\vec{u})c^*~\forall (x,y)\in\Omega~\forall t\in]0;T]\) with the boundary condition \(\vec{u}c\cdot\vec{n}=0~\forall (x,y)\in\partial\Omega\cap \{(x,y)|\vec{u}(x,y,t)\cdot\vec{n}>0\}~\forall t\in]0;T]\).

Type:

AdvectionAdjoint

D

The diffusion linear operator with its adjoint operator \(D^*:c\mapsto \nabla\cdot(\mathbf{K}^T\nabla c)~\forall (x,y)\in\Omega~\forall t\in]0;T]\)

Type:

Diffusion

R

The reaction linear operator with its adjoint operator \(R^*:c\mapsto \tau_{loss}c~\forall (x,y)\in\Omega~\forall t\in]0;T]\).

Type:

Reaction

ID

The identity operator \(Id:c\mapsto c~\forall (x,y)\in\Omega~\forall t\in]0;T]\).

Type:

Id

implemented_solver_type

The list of keywords of the implemented types of time discretization. The implemented types of time discretization are: [‘implicit’, ‘semi-implicit’, ‘semi-implicit with matrix inversion’, ‘implicit with stationnary matrix inversion’]

Type:

list of str

time_discretization

The keyword identifying the type of time discretization. The type of time discretization should be the same the one used for the pheromone propagation model implemented in the class DiffusionConvectionReaction2DEquation to ensure that the adjoint of the discretization of the pheromone propagation model coincides with the discretization of the adjoint model.

Type:

str

tol_inversion

Tolerance for the inversion estimation algorithm called at each time step.

Type:

float

transpose_inv_matrix_implicit_part

Inverse matrix of the implicit part matrix of the time discretization of the adjoint model. If the time discretizations are the same, this matrix coincides with the transpose of the inverse matrix of the implicit part matrix of the time discretization of the pheromone propagation implemented in the class DiffusionConvectionReaction2DEquation. Initialized to None.

Type:

ndarray or None

Notes

In a future version of the module, the adjoint operator of the convection operator \(A^*:c^*\mapsto - \nabla\cdot(\vec{u}c^*)+(\nabla.\vec{u})c^*\) will be implemented in the Advection class, instead of implemented in the AdvectionAdjoint class.

__init__(U, K, coeff_depot, msh, time_discretization='semi-implicit', tol_inversion=1e-14)[source]

Constructor method

Parameters:
  • msh (MeshRect2D) – The geometry of the domain.

  • U (Velocity) – The wind field \(\vec{u}(x,y,t)\).

  • K (DiffusionTensor) – The diffusion tensor \(\mathbf{K}(x,y,t)\).

  • coeff_depot (ndarray) – The deposition coefficient \(\tau_{loss}(x,y)\).

  • obs (Obs) – Object containing all the features related to the observations and estimation of the observed variables

  • time_discretization (str, default: 'semi-implicit') – The keyword identifying the type of time discretization.

  • tol_inversion (float, optional, default: 1e-14) – Tolerance for the inversion estimation algorithm called at each time step.

Raises:

ValueError – if the type of discretization is not implemented, i.e. if time_discretization is not in implemented_solver_type

at_current_time(tc)[source]

Update the attributes A_adjoint and D at a given time.

Parameters:

tc (float or integer) – The current time.

Notes

Updates the attributes A_adjoint and D their own attributes using the method at_current_time() of resp. the class AdvectionAdjoint and the class Diffusion.

init_inverse_matrix(path_to_matrix=None, matrix_file_name=None)[source]

Initialize the attribute transpose_inv_matrix_implicit_part if needed. If the inverse matrix of the implicit part matrix of the time discretization of the associated pheromone propagation model has already been computed, it loads the previously computed inverse matrix and compute the adjoint. Otherwise, raises an error as the matrix should have been computed when initializing the pheromone propagation model as object of the class DiffusionConvectionReaction2DEquation using the method init_inverse_matrix().

Parameters:
  • path_to_matrix (str, optional, default: None) – Path where to save or load the inverse matrix. If not provided, set to ‘./data’.

  • matrix_file_name (str, optional, default: None) – Name of the file where the matrix stored or will be saved. If not provided, set to ‘inv_matrix_**_scheme’ with ** either ‘implicit’ or ‘semi_implicit’ depending on the time discretization.

Raises:

ValueError – if no files correspond to the provided path and file name.

Notes

This method is usefull only if time_discretization is either ‘semi-implicit with matrix inversion’ or ‘implicit with stationnary matrix inversion’. Otherwise, the linear system to solve at each time steps is solved using conjugate gradient or GMRES algorithm, and the attribute transpose_inv_matrix_implicit_part is not used.

onesensor_solver(cost, index_sensor, display_flag=True)[source]

Compute the one-sensor adjoint state \(c^*(x,y,t)\) by solving the adjoint model on the whole time window.

Parameters:
  • cost (Cost) – The cost function, containing especially the gap between the observation data and the prediction of the observation variable \(2\mathbf{R}^{-1}\left(m(c(s))-m_{obs}\right)\).

  • display_flag (bool, default: True) – If True, print the evolution of the solver through the time iterations.

Returns:

p_out – The one-sensor adjoint state \(c^*(x,y,t)\).

Return type:

ndarray

Notes

The output \(c^*(x,y,t)\) has been raveled into a (msh.t_array.size * msh.y.size * msh.x.size,)-shape array to match the format of the attribute value of the class Control.

solver(cost, display_flag=True)[source]

Compute the adjoint state \(c^*(x,y,t)\) by solving the adjoint model on the whole time window.

Parameters:
  • cost (Cost) – The cost function, containing especially the gap between the observation data and the prediction of the observation variable \(2\mathbf{R}^{-1}\left(m(c(s))-m_{obs}\right)\).

  • display_flag (bool, default: True) – If True, print the evolution of the solver through the time iterations.

Returns:

p_out – The adjoint state \(c^*(x,y,t)\).

Return type:

ndarray

Notes

The output \(c^*(x,y,t)\) has been raveled into a (msh.t_array.size * msh.y.size * msh.x.size,)-shape array to match the format of the attribute value of the class Control.

solver_given_adj_deriv_obs_op(adjoint_derivative_obs_operator, cost, display_flag=True)[source]

Compute the adjoint state \(c^*(x,y,t)\) by solving the adjoint model on the whole time window for given the adjoint of the derivative of the observation operator.

Parameters:
  • adjoint_derivative_obs_operator (callable) – Function that returns the adjoint of the derivative of the observation operator map \(\left(\left(\frac{dm}{dc}\right)^*\cdot\delta m\right) (x,y)\) given the current time \(t\) and an element of the space of the observed variable \(\delta m\).

  • cost (Cost) – The cost function, containing especially the gap between the observation data and the prediction of the observation variable \(2\mathbf{R}^{-1}\left(m(c(s))-m_{obs}\right)\).

  • display_flag (bool, default: True) – If True, print the evolution of the solver through the time iterations.

Returns:

p_out – The adjoint state \(c^*(x,y,t)\).

Return type:

ndarray

Notes

The output \(c^*(x,y,t)\) has been raveled into a (msh.t_array.size * msh.y.size * msh.x.size,)-shape array to match the format of the attribute value of the class Control.