InitDAE for computing consistent initial values and integration of DAEs

DAEinitialization module

@author: estevez schwarz/lamour

Copyright (C) 2020

This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see <http://www.gnu.org/licenses/>.

InitDAE.DAEinitialization.consistentValues.JacobianDerivativeArray(X, taylorObject, ti, K, dim)

JacobianDerivativeArray creates the Jacobian matrix of the derivativeArray as two dimensional ndarray

InitDAE.DAEinitialization.consistentValues.consistentValues(example, TaylorOrder, alpha, x0, ti, augmentation=False, ftol=1e-08, disp=True, maxiter=100, info_output='01')
consistentValues computes consistent values of the DAE example in the timepoint ti such that || P(x(ti)-alpha) ||_2 -> min! using the modul minimize with the method SLSQP
param:
  • example - object of example class

  • TaylorOrder - integer, order of the highest derivative of Taylor coefficients

  • alpha - reference vector as ndarray of size dim = dimension of the DAE

  • x0 - ndarray containing initial guess of length dim*TaylorOrder

  • ti - timepoint

  • augmentation
    • if False: using the TaylorOrder
    • if True: augmentation of number of Taylor coefficients starting with 1 up to TaylorOrder+1
  • ftol = 1e-8 - (parameter for python module minimize)

  • disp = True - (parameter for python module minimize)

  • maxiter = 100 - (parameter for python module minimize): maximal number of iterations

  • info_level = 0 - integer, higher level implies more information, 0 … 3

    If info_level = 0 no output with the exception of error messages.

    If info_level > 0 then two log-files are generated: example_name_array_date.log and example_name_consi_date.log

    array contains Jacobian matrices, projectors, ranks, SVD values etc.

    consi contains values related to the computation of consistent initial values and the solution in full accuracy (16 digits)

raised:
  • exceptions:
    • MinimizationNotSuccessful
    • ValueError
return:
  • solution as ndarray of length dim*(TaylorOrder+1) containing values for x(t0), x’(t0),x”(t0)/2!, …, x (TaylorOrder)(t0)/(TaylorOrder)!
    Caution: Only x(t0), x’(t0),x”(t0)/2!, …, x (s)(t0)/s! are consistent for s = TaylorOrder+1-index (see [1])
  • TaylorEx Object of class TaylorDAE of the example (to be used in index_Condition)
InitDAE.DAEinitialization.consistentValues.consistentValuesIndex(example, TaylorOrder, alpha, x0, t0, augmentation=False, ftol=1e-08, disp=False, maxiter=100, info_output='01')
consistentValuesIndex wraps the moduls consistentValues and index_Condition.
param:
  • example - object of example class

  • TaylorOrder - integer, order of the highest derivative of Taylor coefficients

  • alpha - reference vector as ndarray of size dim = dimension of the DAE

  • x0 - ndarray containing initial guess of length dim*TaylorOrder

  • t0 - timepoint

  • augmentation = False
    • if False: using the TaylorOrder
    • if True: augmentation of number of Taylor coefficients starting with 1 up to TaylorOrder+1
  • ftol = 1e-8 - (parameter for python module minimize)

  • disp = False - (parameter for python module minimize)

    If disp = True also results of the iteration of the minimizer are displayed.

  • maxiter = 100 - (parameter for python module minimize): maximal number of iterations

  • info_output - string ‘ab’

  • info_output = higher level implies more information, ‘0’ … ‘3’

    If info_output = ‘0’ no output with the exception of error messages.

    If ‘a’ > ‘0’ then log-file example_name_array_date.log is generated.

    array contains Jacobian matrices, projectors, ranks, SVD values etc.

    If ‘b’ > ‘0’ then log-file example_name_consi_date.log

    consi contains values related to the computation of consistent initial values and the solution in full accuracy (16 digits)

return:
  • solution as ndarray of length dim*(TaylorOrder+1) containing values for x(t0), x’(t0),x”(t0)/2!, …, x (TaylorOrder)(t0)/(TaylorOrder)!
    Caution: Only x(t0), x’(t0),x”(t0)/2!, …, x (s)(t0)/s! are consistent for s = TaylorOrder+1-index (see [1])
  • index of the DAE in the computed consistent point
  • cond - condition number of the index computation
  • computation_robust - If True -> The relation TaylorOrder+1-index-s==0 is fulfilled
InitDAE.DAEinitialization.consistentValues.constraint(taylorObject, ti, K, dim)

constraint delivers the dictionary containing the derivative array of order K of the example (in taylorObject) at timepoint ti.

InitDAE.DAEinitialization.consistentValues.derivativeArray(X, taylorObject, ti, K, dim)

derivativeArray creates the derivativeArray as ndarray

InitDAE.DAEinitialization.consistentValues.func(X, P, alpha, ex)

Objective function ||P(x_0-alpha)||_2, where X = (x_0,x_1,…) :param X - ndarray of Taylor variable :param P - orthogonal projector :param alpha - ‘initial value’ :param ex

InitDAE.DAEinitialization.consistentValues.func_deriv(X, P, alpha, ex)

Derivative of objective function

InitDAE.DAEinitialization.consistentValues.index_Condition(TaylorEx, dimOfDAE, TaylorVector, ti)
index_Condition calculates the index of the DAE and the matrix condition number of the “Jacobian” of the derivative array
param:
  • TaylorEx - object of class TaylorDAE of the example
  • dimOfDAE - dimension of the DAE
  • TaylorVector - TaylorVector as ndarray (e.g. solution of consistentValues)
  • ti - timepoint
return:
  • index - index of the DAE
  • cond - condition number
  • s - level of s-fullness
InitDAE.DAEinitialization.consistentValues.projectorPi(G, P, n)
projectorPi computes the projector Pi - see paper
G - coefficent matrix of constraints P - projector P n - dimension of DAE

Integration module

Created at 20.05.2014/ last change 25.03.2020

@author: estevez/lamour

Copyright (C) 2020

This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see <http://www.gnu.org/licenses/>.

class IntegrationDAE.DAE_TaylorIntegration.TaylorIntegration(example, K, tol=None, info_output='001', method=0)

Bases: InitDAE.DAEanalysis.DerivativeArray.TaylorDAE

From the TaylorDAE class, the TaylorIntegration object inherits the attributes and has the additionally attributes
param:
  • info_output - additional component ‘c’, string ‘abc’ controls the information output, 0 no output up to 3 most information

    a - level for linear algebra, like SVD, rank determination etc.

    b - level for calculation of consistent initial values

    c - level for integration

    standard - ‘001’

  • method - integration method

    = 0: Explicit Taylor K-index=ke, ki=0

    = 1: Fully implicit Taylor K-index=ki, ke=0

    = 2: Two-halfstep-method K-index=ke=ki

    = 3: HOP method for ki=ke

    = 4: HOP method for ki-1=ken

    = 5: HOP method for ki-2=ke

  • FOR INTERNAL USE ONLY

  • Kc - initial value = None

  • ki - initial value = None

  • ke - initial value = None

  • h - initial value = None

  • Xlast - initial value = None

  • HOPweightsE - initial value =None

  • HOPweightsI - initial value =None

  • index - initial value =NaN

HOP_Evaluation(x0, E_I)
Computes xx according to HOP-formulas (e/i)
param:
  • x0 - ndarray

  • h - step size

  • E_I - decides whether we compute x_e or x_i

    for E_I==’E’: x_e

    ‘I’: x_i

  • xx (for x_e or x_i)

Taylor_Method_alpha(Xj, n, h)
Taylor_Method_alpha computes alpha according to the chosen integration method.
param:
  • X - ndarray of Taylor variable
  • n - dimension of DAE
  • h - current stepsize
return:
  • alpha -
Taylor_Method_x0i(X, n)
Compute x0i according to the chosen Taylor Method
param:
  • X - ndarray of Taylor variable
  • n - dimension of DAE
return:
  • x0i
Taylor_Step(x0, h)
Taylor_Step computes an approximation for the Taylor series X0 at t_0+h if the Taylor series X0 at t_0 and h are given.
param:
  • x0(t) - ndarray with Taylor coefficents
  • h - step size
return:
  • y - x0(t+h)
Taylor_integration_adaptiveStep(t0, T, alpha, max_iterations=20, tol=1e-06, hmin=1e-08, K_start=2, disp=False, RobustIntegration=True, ftol=1e-08, hstart=-1.0)

Taylor_integration with adaptive stepsize. The idea consists in the computation of consistent initial values,

the estimation of a relevant stepsize, a step using the Taylor expansion and so on.
param:
  • t0 - initial time point
  • T - end of the integration interval
  • x0 - initial guess of x(t0), see also computation of consistent initial values
  • max_iterations - maximal number of iterations in minimizer
  • tol - rel. integration accuracy
  • hmin - minimal stepsize
  • ftol - accuracy of consistent initial value computation
  • hstart - initial stepsize, if < 0 estimated stepsize is used.
return:
  • xout - ndarray with solution at time points
  • tlist - list of timepoints
  • index - list of DAE index at timepoints
  • condM - list of condition number at timepoints
Taylor_integration_constant_h(t0, x0, T, NumSteps, ftol=1e-14, max_iterations=20)
Taylor__integration_constant_h integrates a DAE over the interval [t0,T] with NumSteps ( = number of intervals )
param:
  • t0 -
  • x0 - initial value
  • T -
  • NumSteps
  • ftol
  • max_iterations
return:
  • xout - ndarray with solution at time points
  • ti_list - list of timepoints
  • index - list of DAE index at timepoints
  • condM - list of condition number at timepoints
Taylor_integration_timeList(ti_list, x0, ftol=1e-14, max_iterations=20)
Taylor__integration_timeList integrates a DAE with given steps in the parameter ti_list
param:
  • ti_list - numpy.ndarray containing [t0,t1,…,ti,…,T]
  • x0 - initial guess
  • ftol = 1e-14
  • max_iterations = 20
return:
  • ti_list - list of successful t points (mostly identical with input list)
  • xout - solution at timepoints in t-_list
  • index - array of DAE index at every timepoint
  • condM - array of condition number at every timepoint
check_solution_along_path(path, ftol=1e-10)
check_solution_along_path delivers the index and the condition number of a DAE along a by a list given path
param:
  • path - String containing the path of a file with list of solutions
    structure of every row of the ASCII-file: timepoint x_1 x_2 … x_n
return:
  • t_list - list of timepoints from file
  • index - list of DAE index at timepoints
  • condition - list of condition numbers at timepoints
check_solutionpath(x, a, b, time_points=100)
check_solutionpath computes the condition of the Taylor coefficient matrix in every point of the solution x in the interval [a,b]
param:
  • x - function of t, which delivers Taylorcoefficents too. (X=algopy.zeros(dimension,dtype=t) )
  • a,b - time interval of interest
  • time_points - number of sampling points
return:
  • condition - ndarray(time_points+1,2), condition[:,0] = list of timepoints, condition[:,1] = list of condition numbers
  • index - ndarray(time_points+1,1) = list of DAE index
consistentStep(TaylorOrder, alpha, x0, ti, augmentation=False, ftol=1e-08, disp=True, maxiter=100, info_output='001')
consistentStep computes consistent values of the DAE example in the timepoint ti such that || P(x(ti)-alpha) ||_2 -> min! using the modul minimize with the method SLSQP
param:
  • self - Taylor Object of the type TaylorIntegration

  • TaylorOrder - integer, order of the highest derivative of Taylor coefficients

  • alpha - reference vector as ndarray of size dim = dimension of the DAE

  • x0 - ndarray containing initial guess of length dim*TaylorOrder

  • ti - timepoint

  • augmentation
    • if False: using the TaylorOrder
    • if True: augmentation of number of Taylor coefficients starting with 1 up to TaylorOrder+1
  • ftol = 1e-8 - (parameter for python module minimize)

  • disp = True - (parameter for python module minimize)

  • maxiter = 100 - (parameter for python module minimize): maximal number of iterations

  • info_output - additional component ‘c’, string ‘abc’ controls the information output, 0 no output up to 3 most information

    a - level for linear algebra, like SVD, rank determination etc.

    b - level for calculation of consistent initial values

    c - level for integration

  • info_level = 0 - integer, higher level implies more information, 0 … 3

    If info_level = 0 no output with the exception of error messages.

    If info_level > 0 then two log-files are generated: example_name_array_date.log and example_name_consi_date.log

    array contains Jacobian matrices, projectors, ranks, SVD values etc.

    consi contains values related to the computation of consistent initial values and the solution in full accuracy (16 digits)

raised:
  • exceptions:
    • MinimizationNotSuccessful
    • ValueError
return:
  • solution as ndarray of length dim*(TaylorOrder+1) containing values for x(t0), x’(t0),x”(t0)/2!, …, x (TaylorOrder)(t0)/(TaylorOrder)!
    Caution: Only x(t0), x’(t0),x”(t0)/2!, …, x (s)(t0)/s! are consistent for s = TaylorOrder+1-index (see [1])
  • TaylorEx Object of class TaylorDAE of the example (to be used in index_Condition)
consistentStepIndex(TaylorOrder, alpha, x0, t0, augmentation=False, ftol=1e-08, disp=False, maxiter=100, info_output='001')
consistentStepIndex wraps the moduls consistentValues and index_Condition.
param:
  • self - Taylor Object of the type TaylorIntegration

  • TaylorOrder - integer, order of the highest derivative of Taylor coefficients

  • alpha - reference vector as ndarray of size dim = dimension of the DAE

  • x0 - ndarray containing initial guess of length dim*TaylorOrder

  • t0 - timepoint

  • augmentation = False
    • if False: using the TaylorOrder
    • if True: augmentation of number of Taylor coefficients starting with 1 up to TaylorOrder+1
  • ftol = 1e-8 - (parameter for python module minimize)

  • disp = False - (parameter for python module minimize)

    If disp = True also results of the iteration of the minimizer are displayed.

  • maxiter = 100 - (parameter for python module minimize): maximal number of iterations

  • info_output - additional component ‘c’, string ‘abc’ controls the information output, 0 no output up to 3 most information

    a - level for linear algebra, like SVD, rank determination etc.

    b - level for calculation of consistent initial values

    c - level for integration

return:
  • solution as ndarray of length dim*(TaylorOrder+1) containing values for x(t0), x’(t0),x”(t0)/2!, …, x (TaylorOrder)(t0)/(TaylorOrder)!
    Caution: Only x(t0), x’(t0),x”(t0)/2!, …, x (s)(t0)/s! are consistent for s = TaylorOrder+1-index (see [1])
  • index of the DAE in the computed consistent point
  • cond - condition number of the index computation
  • computation_robust - If True -> The relation TaylorOrder+1-index-s==0 is fulfilled
static func_int(X, P, alpha, Object)
Objective function ||P(x_0-alpha)||_2, where X = (x_0,x_1,…)
param:
  • X - ndarray of Taylor variable
  • P - orthogonal projector
  • alpha - initial guess
  • Object - Taylor object
return:
  • sqrt(yty) - ||y||_2 with y = P(x0-alpha)
static func_int_deriv(X, P, alpha, Object)
Derivative of objective function
param:
  • X - ndarray of Taylor variable
  • P - orthogonal projector
  • alpha - ‘initial value’
  • Object - Taylor object
return:
  • deriv -
logging_integration_solution(t, x)

logging_integration_solution writes a tuple of t,x(t) to a logging file

plot_figure(tlist, x, mylabel, xlab='t', ylog=False)
plot_figure supports the output of graphics of results stored in ndarrays of the form [ number of timepoints, dim of the DAE]
param:
  • tlist - list of timepoints
  • x - list of solution at the timepoints
  • mylabel - string for labeling x
  • xlab - label of x-axis
  • ylog - switch for logarithmic abscissa
return:

no output, plot opens

weights_HOP()
IntegrationDAE.DAE_TaylorIntegration.adapt_solution_list(path, K)

adapt_solution_list computes a list of Taylor function values up to order K, where the values of the

solution in the file “path” are used and the first derivative approximated by finite differences.

The other Taylor coefficients from 2 to K are set to zero.
param:
  • path - (relative) path to the file with solution in the form t, x(t)
  • K - number of Taylor coefficients
return:
  • Taylor object with solution, used in check_solution_along_path
IntegrationDAE.DAE_TaylorIntegration.compute_finite_differences(solution_list)

compute_finite_differences delivers a list of t_i, x(t_i), x’(t_i).

The first derivative is approximated by finite differences of second order
param:
  • solution_list - array with columns of time point, x(1),…x(M-1)
return:
  • derivative_list - ndarray with x(t_i), x’(t_i), used in adapt_solution_list

DAEanalysis module

Created on 20.05.2014/ last change 23.11.19

@author: estevez schwarz/lamour

Copyright (C) 2020

This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see <http://www.gnu.org/licenses/>.

class InitDAE.DAEanalysis.DerivativeArray.TaylorDAE(example, K, tol=None, info_output='01')
TaylorDAE object has the attributes
attribute:
  • example - object of type “example” (see Implementation of an example)

  • K - Number of Taylor coefficients 0,…,K, i.e., we have K+1 coefficients

  • tol - absolute accuracy for rank determination

  • info_output - string ‘ab’ controls the information output, 0 no output up to 3 most information

    a - level for linear algebra, like SVD, rank determination etc.

    b - level for calculation of consistent initial values

    standard - ‘01’

CheckOmegaTaylorAndProjectors(Coeff_Matrix, n)
CheckOmegaTaylorAndProjectors checks the index of the DAE and computes the projectors T_j.
param:
  • Coeff_Matrix - Jacobian matrix of the derivative array with additional A at the first block row
  • n - Number of equations of the original DAE
return:
  • T - Projectors T_j, all in one matrix
  • index - Differentiation index of the DAE
Check_s_full(CoeffMatrix, n)

Check_s_full computes the value s, which characterizes the accuracy of the result.

IndexAndCondition(tti, zt)
IndexAndCondition computes the index and the condition number of the DAE.
param:
  • tti - Taylor object of ti
  • zt - Taylor object of z(t)
return:
  • index - index of the DAE
  • condM - condition number of index determination
  • s - level of s-fullness
Jacobian_forCoefficentMatrix(tti, zt)
The method Jacobian_forCoefficentMatrix computes the Jacobian matrices Atilde and Btilde.
param:
  • tti - Taylor object of timepoint
  • zt - MyUTPM object of z
return:
  • The output depends of the structure of the DAE. For
    • type_of_DAE = ‘proper’
      • Atilde= A.D, Btilde= A D’+ B with f( d’(x,t),x,t) and
        A = df/dx’(zt,ti), D = dd/dx(zt,ti), B = df/dx(zt,ti)
    • type_of_DAE = ‘standard’
      • Atilde= A, Btilde= B with f( x’,x,t)
    • type_of_DAE = ‘linear’
      • Atilde= eval_A(t), Btilde= eval_B(t) with f( x’,x,t) = eval_A x’+eval_B x-q
raised:
  • Exception - unknown type of DAE
Jacobian_matrices_AB(ti, zt)
The method Jacobian_matrices_AB computes the Jacobian matrices
param:
  • ti - Taylor object
  • zt - MyUTPM object
return:
  • A = df/dx’(zt,ti)
  • B = df/dx(zt,ti)
  • f_Taylor - derivative array of f (see MyUTPM.tracing_A_B)
Jacobian_matrices_ADB(ti, zt)
The method Jacobian_matrices_ADB computes the Jacobian matrices
param:
  • ti - Taylor object
  • zt - MyUTPM object
return:
  • A = df/dx’(zt,ti)
  • D = dd/dx(zt,ti)
  • B = df/dx(zt,ti)
Jacobian_matrix_u(zt)
The method Jacobian_matrix_u computes the Jacobian matrix of u(x)
param:
  • zt - MyUTPM object
return:
  • u_x = du/dx(zt)
coefficent_matrix(tti, xt)
coefficient matrix of Taylor method
param:
  • tti - Taylor representation of ti (time point)
  • xt - Taylor coefficients of x(ti)
return:
  • Coeff_Matrix - see structure below
  • block_dimension - n
structure of Coeff_Matrix for constant standard DAE Ax’+Bx = q
A          
B A       n(K+1) x n(K+1)
  B A     self.K+1 - number of Taylor coefficients
      A, B nxn matrices
      B A  
coefficent_matrix_optimization(tti, xt)
coefficient matrix of Taylor method
param:
  • tti - Taylor representation of ti (time point)
  • xt - Taylor coefficients of x(ti)
return:
  • Coeff_Matrix -
  • block_dimension -
structure of Coeff_Matrix for constant standard DAE Ax’+Bx = q
B A       n*K x n(K+1)
  B A     self.K+1 - number of Taylor coefficients
      A, B nxn matrices
      B A  
f_taylor(x, t)
f_taylor computes the Taylor coefficients of f depending of the DAE type
param:
  • x - Taylor object of x(t)
  • t - Taylor object of t
return:
  • f - Taylor expansion of f
raised:
  • Exception - unknown type of DAE
log_files()

log_files provides a list which contains the Logger for the different files.

The following numerical submethods are logged:
  • 0 - computation of consistent initial values
  • 1 - Newton method (not active for InitDAE)
  • 2 - linear algebra results of derivative array computations
  • 3 - integration of IVPs (not active for InitDAE)
static write_Values2file(name, message, dim, r, S1, Sr, Sr1, Send, tol, info_level)
write_Values2file writes the values of the rank determination to a file
param:
  • name
  • message
  • dim
  • r - see toolboxLA.svdinfo.CompCondSing
  • S1
  • Sr
  • Sr1
  • Send
  • tol
  • info_level - if >‘2’ information printed on terminal

Created on 13.01.2016

@author: lamour

exception InitDAE.DAEanalysis.NumericalException.FNotComputable

Bases: Exception

exception InitDAE.DAEanalysis.NumericalException.IncreasingTaylorCoefficient(K_start, K, K_end)

Bases: Exception

exception InitDAE.DAEanalysis.NumericalException.LittleNumberOfTaylorCoefficients(K_start)

Bases: Exception

exception InitDAE.DAEanalysis.NumericalException.MinimalDampingParameter

Bases: Exception

exception InitDAE.DAEanalysis.NumericalException.MinimalStepsize

Bases: Exception

exception InitDAE.DAEanalysis.NumericalException.MinimizationNotSuccessful

Bases: Exception

exception InitDAE.DAEanalysis.NumericalException.SingularDAE(string)

Bases: Exception

toolboxLA module

Created on Thu Mar 05 09:07:30 2015

@author: estevez schwarz

Copyright (C) 2020

This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see <http://www.gnu.org/licenses/>.

InitDAE.toolboxLA.svdinfo.CompCondSing(M, tol=None)
CompCondSing computes the ratio between the largest and the smallest nonzero singular value of a matrix M and returns some related information.
param:
  • M - Matrix M
  • tol - tolerance for the rank decision (optional). The computation of the default tolerance corresponds to linalg.matrix_rank
return:
  • cond - S[0]/S[rankM-1]: ratio between the largest and the smallest nonzero singular value of M
  • rankM - considered rank of the matrix M
  • S[0] - largest singular value
  • S[rankM-1] - smallest singular value considered as positive
  • S[min(len(S)-1,rankM)] - largest singular value considered as zero
  • S[-1] - smallest singular value
  • tol - used tolerance
InitDAE.toolboxLA.svdinfo.OrthogonalBasisProjectorAlongIm(A, tol=None)
OrthogonalBasisProjectorAlongIm computes the orthogonal basis to build the projector along im(A) for a given matrix A
param:
  • A - Matrix A
return:
  • Ur.T - resulting orthogonal basis
  • rS - considered rank of the matrix A (and corresponding S from SVD)
  • S[0] - largest singular value of A
  • S[rS-1] - smallest singular value considered as positive
  • S[min(len(S)-1,rS)] - largest singular value considered as zero
  • S[-1] - smallest singular value
InitDAE.toolboxLA.svdinfo.OrthogonalBasisProjectorOntoKer(A, tol=None)
OrthogonalBasisProjectorOntoKer computes the orthogonal basis to build the projector onto ker(A) for a given matrix A
param:
  • A - Matrix A
return:
  • VTr - resulting orthogonal basis
  • rS - considered rank of the matrix A (and corresponding S from SVD)
  • S[0] - largest singular value of A
  • S[rS-1] - smallest singular value cosidered as positive
  • S[min(len(S)-1,rS)] - largest singular value considered as zero
  • S[-1] - smallest singular value
InitDAE.toolboxLA.svdinfo.OrthogonalProjectorAlongIm(A, tol=None)
OrthogonalProjectorAlongIm computes the orthogonal projector along im(A) for a given matrix A
param:
  • A - Matrix A
return:
  • W - resulting orthogonal projector
  • rS - considered rank of the matrix A (and corresponding S from SVD)
  • S1 - largest singular value of A
  • Sr - smallest singular value considered as positive
  • Sr1 - largest singular value considered as zero
  • Send - smallest singular value
InitDAE.toolboxLA.svdinfo.OrthogonalProjectorOntoKer(A, tol=None)
OrthogonalProjectorOntoKer computes the orthogonal projector onto ker(A) for a given matrix A
param:
  • A - Matrix A
return:
  • Q - resulting orthogonal projector
  • rS - considered rank of the matrix A (and corresponding S from SVD)
  • S1 - largest singular value of A
  • Sr - smallest singular value cosidered as positive
  • Sr1 - largest singular value considered as zero
  • Send - smallest singular value
InitDAE.toolboxLA.svdinfo.RankSigma(S, tol=None)
RankSigma computes the rank of the diagonal matrix Sigma (S) from the SVD
param:
  • S - Matrix S
  • tol - tolerance for the rank decision (optional). The computation of the default tolerance corresponds to linalg.matrix_rank
return:
  • rank - rank(S)
  • tol - used tolerance

toolboxAD module

Created on 13.01.2016

@author: lamour

Copyright (C) 2020

This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see <http://www.gnu.org/licenses/>.

class InitDAE.toolboxAD.MyUTPM.MyUTPM(X)

Bases: algopy.utpm.utpm.UTPM

A_B(cgAB, f, t)
Calculation of the Jacobians A and B
param:
  • cgAB - computational graph
  • f - function f
  • t - UTPM object - taylor_init of t
return:
  • A - Taylor object of Jacobian df/dx’
  • B - Taylor object of Jacobian df/dx
A_D_B_new(cg, cgAB, d, f, t)
Calculation of the Jacobians A, D and B
param:
  • cg - computational graph of f
  • cgAB - computational graph of A, B
  • d - function d of proper formulation
  • f - function f
  • t - UTPM object - taylor_init of t
return:
  • D - Taylor object of Jacobian matrix dd/dx
  • A - Taylor object of Jacobian matrix df(z,x,t)/dz
  • B - Taylor object of Jacobian matrix df(z,x,t)/dx
diff()

differentiation for UTPM objects by shifting the coefficients

static taylor_init(D, t0)

initialization of independent variable t at t0

tracing_A_B(f, t)
The function traces the function f to compute the Jacobians A and B
param:
  • f - function f
  • t - time point
return:
  • cgAB - computational graph of f
  • Fyxtt.x.data - ndarray of Taylorcoefficients of f
tracing_A_D_B(d, f, t)

The function traces the functions d and f to compute the Jacobians D, A and B