InitDAE for computing consistent initial values¶
DAEinitialization module¶
@author: estevez schwarz/lamour
Copyright (C) 2018
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/>.
-
DAEinitialization.consistentValues.
JacobianDerivativeArray
(X, taylorObject, ti, K, dim)¶ JacobianDerivativeArray creates the Jacobian matrix of the derivativeArray as two dimensional ndarray
-
DAEinitialization.consistentValues.
consistentValues
(example, TaylorOrder, alpha, x0, ti, augmentation=False, ftol=1e-08, disp=True, maxiter=100, info_level=3)¶ - 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)
-
DAEinitialization.consistentValues.
consistentValuesIndex
(example, TaylorOrder, alpha, x0, t0, augmentation=False, ftol=1e-08, disp=False, maxiter=100, info_level=0)¶ - 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_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)
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
-
DAEinitialization.consistentValues.
constraint
(taylorObject, ti, K, dim)¶ constraint delivers the dictionary containing the derivative array of order K of the example (in taylorObject) at timpoint ti.
-
DAEinitialization.consistentValues.
derivativeArray
(X, taylorObject, ti, K, dim)¶ derivativeArray creates the derivativeArray as ndarray
-
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
-
DAEinitialization.consistentValues.
func_deriv
(X, P, alpha, ex)¶ Derivative of objective function
-
DAEinitialization.consistentValues.
index_Condition
(Example, TaylorEx, dimOfDAE, TaylorVector, ti, info_level=0)¶ - index_Condition calculates the index of the DAE and the matrix condition of the “Jacobian” of the derivative array
param: - Example - object of the example
- 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
DAEanalysis module¶
Created on 20.05.2014/ last change 29.1.18
@author: estevez schwarz/lamour
-
class
DAEanalysis.DerivativeArray.
TaylorDAE
(example, K, tol=None, info_level=2)¶ - 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_level controls the information output
- 0 - no output
- 1 ... 4 output increases
-
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
-
static
write_Values2logger
(logger, info_level, print_level, message)¶ write_Values2logger stores messages depending of the info_level in the logger.
param: - logger - used logger
- info_level - if > 0 creates an output. The greater the more.
- print_level - determines the level of the message.
- message - string containing the message
Created on 13.01.2016
@author: lamour
-
exception
DAEanalysis.NumericalException.
FNotComputable
¶ Bases:
exceptions.Exception
-
exception
DAEanalysis.NumericalException.
IncreasingTaylorCoefficient
(K_start, K, K_end)¶ Bases:
exceptions.Exception
-
exception
DAEanalysis.NumericalException.
LittleNumberOfTaylorCoefficients
(K_start)¶ Bases:
exceptions.Exception
-
exception
DAEanalysis.NumericalException.
MinimalDampingParameter
¶ Bases:
exceptions.Exception
-
exception
DAEanalysis.NumericalException.
MinimizationNotSuccessful
¶ Bases:
exceptions.Exception
-
exception
DAEanalysis.NumericalException.
SingularDAE
(string)¶ Bases:
exceptions.Exception
toolboxLA module¶
Created on Thu Mar 05 09:07:30 2015
@author: estevez schwarz
-
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
-
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
-
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
-
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
-
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
-
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
-
class
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
-