Reference

Contents

Index

Base.:==Method
(==)(ps1::T, ps2::T) where T <: AbstractPartitionedStructure

Return true if both partitioned-structures are composed of the same amont of element-structures, and have the same size.

source
PartitionedStructures.M_abstract_part_struct.check_epv_epmMethod
bool = check_epv_epm(epm::Y, epv::Z) where {Y <: AbstractPartitionedStructure, Z <: AbstractPartitionedStructure}

Similar to ==, but it can compare different partitioned-structures, example: an Elemental_pv and an Elemental_pm. check_epv_epm is a superficial test, see full_check_epv_epm(epm, epv) for a complete check of the partitioned-structure (i.e. if each element depends of the same variable subset).

source
PartitionedStructures.M_abstract_part_struct.get_component_listMethod
list = get_component_list(ps::T) where T <: AbstractPartitionedStructure
ith_component = get_component_list(ps::T, i::Int) where T <: AbstractPartitionedStructure

Return either the list of every element-structure composing the partitioned-structure ps or the i-th element-structure of ps.

source
PartitionedStructures.M_abstract_part_struct.get_ee_structMethod
ee_vector = get_ee_struct(eps::AbstractPartitionedStructure{T}) where T
ee = get_ee_struct(eps::AbstractPartitionedStructure{T}, i::Int) where T

Return a vector of every elemental elements fom elemental partitioned-structure eps or only its i-th elemental element.

source
PartitionedStructures.M_abstract_element_struct.get_indicesMethod
indices = get_indices(elt::T) where T <: Element_struct
indice = get_indices(elt::T, i::Int) where T <: Element_struct

Every element-structure is based on a variable subset of a partitioned-structure. get_indices(elt) retrieves the variable set of an element elt. get_indices(elt, i) retrieves the i-th variable associated to elt.

source
PartitionedStructures.M_elt_vec.set_vec!Method
set_vec!(ev::T, vec::Vector{Y}) where {Y <: Number, T <: Elt_vec{Y}}
set_vec!(ev::T, i::Int, val::Y) where {Y <: Number, T <: Elt_vec{Y}}

Set ev.vec to vec or ev.vec[i] = val of the element-vector ev.

source
PartitionedStructures.ModElemental_pv.create_epvMethod
epv = create_epv(sp_set::Vector{SparseVector{T,Y}}; kwargs...) where {T,Y}
epv = create_epv(eev_set::Vector{Elemental_elt_vec{T}}; n=max_indices(eev_set)) where T

Create an elemental partitioned-vector from a vector eev_set of: SparseVector, elemental element-vector or a vector of indices.

source
PartitionedStructures.ModElemental_pv.epv_from_vMethod
epv = epv_from_v(x::Vector{T}, shape_epv::Elemental_pv{T}) where T

Define a new elemental partitioned-vector from x that have the same structure than shape_epv. The value of each elemental element-vector comes from the corresponding indices of x. Usefull to define Uᵢ x, ∀ x.

source
PartitionedStructures.ModElemental_pv.get_eev_subsetMethod
eev_subset = get_eev_subset(epv::Elemental_pv{T}, indices::Vector{Int}) where T

Return a subset of the elemental element vector composing the elemental partitioned-vector epv. indices selects the differents elemental element-vector needed.

source
PartitionedStructures.ModElemental_pv.get_eev_valueMethod
eev_i_value = get_eev_value(epv::Elemental_pv{T}, i::Int) where T
eev_ij_value = get_eev_value(epv::Elemental_pv{T}, i::Int, j::Int) where T

Return either the value of the i-th elemental element-vector of the elemental partitioned-vector epv or only the j-th component of the i-th elemental element-vector.

source
PartitionedStructures.ModElemental_pv.part_vecMethod
epv = part_vec(;n::Int=9, T=Float64, nie::Int=5, overlapping::Int=1, mul::Float64=1.)

Define an elemental partitioned-vector formed by N (deduced from n and nie) elemental element-vectors of size nie. Each elemental element-vector overlaps the previous and the next element by overlapping.

source
PartitionedStructures.ModElemental_pv.prod_part_vectorsMethod
(acc, res) = prod_part_vectors(epv1::Elemental_pv{T}, epv2::Elemental_pv{T}) where T

Perform an elementwise scalar product between the two elemental partitioned-vector epv1 and epv2. acc accumulates the sum of the element-vectors scalar product. res contrains the details of every element-vector scalar product.

source
PartitionedStructures.ModElemental_pv.scale_epv!Method
v = scale_epv!(epv::Elemental_pv{T}, scalars::Vector{T}) where T

Return a vector v from epv where the contribution of each element-vector is multiply by the corresponding value from scalars. v is extract from epv.v

source
PartitionedStructures.ModElemental_pv.set_eev!Method
set_eev!(epv::Elemental_pv{T}, i::Int, vec::Vector{T}) where T
set_eev!(epv::Elemental_pv{T}, i::Int, j::Int, val:: T) where T

Set either the i-th elemental element-vector epv to vec or its j-th component to val.

source
PartitionedStructures.M_elt_mat.get_current_untouchedMethod
index = get_current_untouched(cem::Counter_elt_mat)

Return index: the number of the last partitioned-updates that did not update the element-matrix elt_mat. If the last partitioned-update updates elt_mat then index will be equal to 0.

source
PartitionedStructures.M_elt_mat.get_indexMethod
index = get_index(elt_mat::T) where T <: Elt_mat

Return index: the number of the last partitioned-updates that did not update the element-matrix elt_mat. If the last partitioned-update updates elt_mat then index will be equal to 0.

source
PartitionedStructures.M_elt_mat.iter_infoMethod
(current_update, current_untouched, current_reset) = iter_info(cem::Counter_elt_mat)

Return the information about the last partitioned quasi-Newton update applied onto the element counter cem (associated to an element-matrix).

source
PartitionedStructures.M_elt_mat.total_infoMethod
(total_update, total_untouched, current_reset) = iter_info(cem::Counter_elt_mat)

Return the informations about all the quasi-Newton updates applied onto the element associated to the counter cem.

source
PartitionedStructures.M_part_mat.get_eelo_setMethod
eelo_set = get_eelo_set(plm::T) where T <: Part_LO_mat
eelo = get_eelo_set(plm::T, i::Int) where T <: Part_LO_mat

Return either the vector of every elemental element linear-operator plm.eelo_set or the i-th elemental element linear-operator plm.eelo_set[i].

source
PartitionedStructures.M_part_mat.get_eelo_sub_setMethod
eelo_subset = get_eelo_sub_set(plm::T, indices::Vector{Int}) where T <: Part_LO_mat

Return a subset of the elemental element linear-operators composing the elemental partitioned limited-memory operator plm. indices selects the differents elemental element linear-operators needed.

source
PartitionedStructures.M_part_mat.get_spmMethod
spm = get_spm(pm::T) where T <: Part_mat
spm_ij = get_spm(pm::T, i::Int, j::Int) where T <: Part_mat

Get either the sparse matrix associated to the partitioned-matrix pm or pm[i,j].

source
PartitionedStructures.M_part_mat.set_spm!Method
set_spm!(pm::P) where {T <: Number, P <: Part_mat{T}}

Build the sparse matrix of the partitioned-matrix pm in pm.spm by gathering the contribution of every element-matrix. The sparse matrix is built with respect to the indices of each elemental element linear-operator.

source
PartitionedStructures.M_part_v.add_v!Method
add_v!(pv::T, i::Int, value::Y) where {Y, T <: Part_v{Y}}
add_v!(pv::T, indices::Vector{Int}, values::Vector{Y}) where {Y, T <: Part_v{Y}}

Add value (resp values) to the vector of the partitioned-vector pv.v at the indice i (resp indices).

source
PartitionedStructures.M_part_v.set_v!Method
set_v!(pv::T, v::Vector{Y}) where {Y, T <: Part_v{Y}}
set_v!(pv::T, i::Int, value::Y) where {Y, T <: Part_v{Y}}

Set the components of the vector pv.v (resp. pv.v[i]) from the partitioned-vector pv to the vector v (resp. value).

source
PartitionedStructures.ModElemental_pm.Elemental_pmType
Elemental_pm{T} <: Part_mat{T}

Represent an elemental partitioned quasi-Newton linear-operator. Each element is an elemental element-matrix which may apply a BFGS or a SR1 update. N is the number of elements. n is the size of the elemental partitioned limited-memory operator. eelo_set is the set of elemental element linear-operators. spm and L are sparse matrices either to form the sparse matrix gathering the elements or the Cholesky factor of spm. component_list summarizes for each variable i (∈ {1,..., n}) the list of elements (⊆ {1,...,N}) being parametrised by i. permutation is the current permutation of the elemental partitioned limited-memory operator ([1:n;] initially).

source
Base.permute!Method
permute!(epm::Elemental_pm{T}, p::Vector{Int}) where T

Apply the permutation p to the elemental partitionned matrix epm. The permutation is applied to every elemental element-matrix eem via indices. The current epm permutation is stored in epm.permutation.

source
PartitionedStructures.ModElemental_pm.identity_epmMethod
epm = identity_epm(N::Int, n::Int; T=Float64, nie::Int=5, linear_vector::Vector{Bool} = zeros(Bool, N))

Return a partitionned matrix of type T of N identity elemental element-matrices. Each elemental element-matrix is of size nie with randoms positions. linear_vector indicates which element linear-opeartor should not contribute to the partitioned linear-operator.

source
PartitionedStructures.ModElemental_pm.identity_epmMethod
epm = identity_epm(element_variables::Vector{Vector{Int}}; N::Int=length(element_variables), n::Int=max_indices(element_variables), T=Float64, linear_vector::Vector{Bool} = zeros(Bool, N))
epm = identity_epm(element_variables::Vector{Vector{Int}}, N::Int, n::Int; T=Float64, linear_vector::Vector{Bool} = zeros(Bool, N))

Return a partitionned matrix of type T of N identity elemental element-matrices. N and n may be extrapolate from element_variables. The elemental variables are based from the indices informed in element_variables. linear_vector indicates which element linear-opeartor should not contribute to the partitioned linear-operator.

source
PartitionedStructures.ModElemental_pm.n_i_SPSMethod
epm = n_i_SPS(n::Int; T=Float64, nie::Int=5, overlapping::Int=1, mul=5.)

Define an elemental partitioned-matrix epm of size n. The partitioned-matrix is composed by N ≈ (n/nie)*2 elemental element-matrices, of size nie, they overlap onto the next element by overlapping. The diagonal terms of each elemental element-matrix are of value mul, whereas the other terms are set to 1.

source
PartitionedStructures.ModElemental_pm.n_i_sepMethod
epm = n_i_sep(n::Int; T=Float64, nie::Int=5, mul=5.)

Define an elemental partitioned-matrix epm composed of nie separable blocs. Each elemental element-matrix is composed of 1 except the diagonal terms which are of value mul.

source
PartitionedStructures.ModElemental_pm.ones_epmMethod
epm = ones_epm(N::Int, n::Int; T=Float64, nie::Int=5)

Create a partitionned matrix of type T of N elemental element-matrices ones(nie, nie) whose positions are random. The partitionned matrix created may be singular.

source
PartitionedStructures.ModElemental_pm.ones_epm_and_idMethod
epm = ones_epm_and_id(N::Int, n::Int; T=Float64, nie::Int=5)

Create a partitionned matrix of type T with N+n elemental element-matrices. The first N elemental element-matrices are ones(nie, nie) with randoms positions in the range 1:n. The remaining n elemental element-matrices are of size 1, with value [1], they are placed in the diagonal terms This way, the partitionned matrix is generally not singular.

source
PartitionedStructures.ModElemental_pm.part_matMethod
epm = part_mat(;n::Int=9, T=Float64, nie::Int=5, overlapping::Int=1, mul=5.)

Define an elemental partitioned-matrix epm composed of N (deduced from n and nie) elemental element-matrices of size nie. Each elemental element-matrix overlaps the previous and the next element by overlapping.

source
PartitionedStructures.ModElemental_em.Elemental_emType
Elemental_em{T} <: DenseEltMat{T}

Represent an elemental element-matrix. It has fields:

  • indices: indices of elemental variables;
  • nie: elemental size (=length(indices));
  • Bie::Symmetric{T, Matrix{T}}: the elemental matrix;
  • counter: counts how many update the elemental matrix goes through from its allocation;
  • convex: if convex==true, then Elemental_em default update is BFGS otherwise it is SR1;
  • linear: if linear==true, then the element matrix contribution is null;
  • _Bsr: a vector used during quasi-Newton update of the elemental matrix.
source
Base.permute!Method
permute!(eem::Elemental_em{T}, p::Vector{Int}) where T

Set the indices of the element variables of eem to p. Must be use with caution.

source
PartitionedStructures.ModElemental_em.fixed_ones_eemMethod
eem = fixed_ones_eem(i::Int, nie::Int; T=Float64, mul=5.)

Create a nie elemental element-matrix of type T at indices index:index+nie-1. All the components of the element-matrix are set to 1 except the diagonal terms that are set to mul. This method is used to define diagonal dominant element-matrix.

source
PartitionedStructures.Utils.BFGS!Method
BFGS!(x0::Vector{Y}, x1::Vector{Y}, g0::Vector{Y}, g1::Vector{Y}, B0::Array{Y,2}; kwargs...) where Y <: Number
BFGS!(s::Vector{Y}, y::Vector{Y}, B::Symmetric{Y,Matrix{Y}}; kwargs...) where Y <: Number
BFGS!(s::Vector{Y}, y::Vector{Y}, B::Array{Y,2}; index=0, reset=4, kwargs...)

Perform the BFGS update in place of the matrix B1 by using the vectors s = x1 - x0 and y = g1 - g0 and the current matrix B0.

source
PartitionedStructures.Utils.BFGSMethod
BFGS(s::Vector{Y}, y::Vector{Y}, B::Array{Y,2}; kwargs...) where Y <: Number
BFGS(x0::Vector{Y}, x1::Vector{Y}, g0::Vector{Y}, g1::Vector{Y}, B::Array{Y,2}; kwargs...) where Y <: Number

Perform the BFGS update over the matrix B by using the vectors s = x1 - x0 and y = g1 - g0.

source
PartitionedStructures.Utils.SE!Method
SE!(s::Vector{Y}, y::Vector{Y}, B::Array{Y,2}; kwargs...) where Y <: Number
SE!(x::Vector{Y}, x1::Vector{Y}, g::Vector{Y}, g1::Vector{Y}, B::Array{Y,2}; kwargs...) where Y <: Number
SE!(s::Vector{Y}, y::Vector{Y}, B::Array{Y,2}, B1::Array{Y,2}; index=0, reset=4, ω = 1e-6, kwargs...)

Perform a BFGS update in place of B1 by using the matrix B, the vectors s = x1 - x0 and y = g1 - g0 if the curvature condition dot(s,y) > eps(eltype(s)) holds. Otherwise, it performs a SR1 update onto B1 with B, s, y.

source
PartitionedStructures.Utils.SEMethod
SE(s::Vector{Y}, y::Vector{Y}, B::Array{Y,2}; kwargs...) where Y <: Number
SE(x::Vector{Y}, x1::Vector{Y}, g::Vector{Y}, g1::Vector{Y}, B::Array{Y,2}; kwargs...) where Y <: Number

Perform a BFGS update over the matrix B by using the vectors s = x1 - x0 and y = g1 - g0 if the curvature condition dot(s,y) > eps(eltype(s)) holds. Otherwise, it performs a SR1 update with B, s, y.

source
PartitionedStructures.Utils.SR1!Method
SR1!(s::Vector{Y}, y::Vector{Y}, B::Array{Y,2}; kwargs...) where Y <: Number
SR1!(x::Vector{Y}, x1::Vector{Y}, g::Vector{Y}, g1::Vector{Y}, B::Array{Y,2}; kwargs...) where Y <: Number
SR1!(s::Vector{Y}, y::Vector{Y}, B::Array{Y,2}, B1::Array{Y,2}; index=0, reset=4, ω = 1e-6, kwargs...)

Perform the SR1 update in place of the matrix B1 by using the vectors s = x1 - x0 and y = g1 - g0 and the current matrix B.

source
PartitionedStructures.Utils.SR1Method
SR1(s::Vector{Y}, y::Vector{Y}, B::Array{Y,2}; kwargs...) where Y <: Number
SR1(x::Vector{Y}, x1::Vector{Y}, g::Vector{Y}, g1::Vector{Y}, B::Array{Y,2}; kwargs...) where Y <: Number

Perform the SR1 update over the matrix B by using the vectors s = x1 - x0 and y = g1 - g0.

source
PartitionedStructures.Utils.max_indicesMethod
indice_max = max_indices(list_of_element_variables::Vector{Vector{T}}) where T
indice_max = max_indices(elt_set::Vector{T}) where T <: Element_struct

Return the maximum index of the element variables in list_of_element_variables or in elt_set.

source
PartitionedStructures.Utils.min_indicesMethod
indice_min = min_indices(list_of_element_variables::Vector{Vector{T}}) where T
indice_min = min_indices(elt_set::Vector{T}) where T <: Element_struct

Return the minimum index of the element variables in list_of_element_variables or in elt_set.

source
PartitionedStructures.Link.eplo_lbfgs_from_epvMethod
eplo = eplo_lbfgs_from_epv(epv::T) where {Y <: Number, T <: Elemental_pv{Y}}

Create an elemental limited-memory partitioned quasi-Newton operator PLBFGS eplo with the same partitioned structure than epv. Each element linear-operator of eplo is set to a LBFGSOperator of suitable size.

source
PartitionedStructures.Link.eplo_lose_from_epvMethod
eplo = eplo_lose_from_epv(epv::Elemental_pv{T}) where {T <: Number}

Create an elemental limited-memory partitioned quasi-Newton operator PLSE eplo with the same partitioned structure than epv. Each element linear-operator of eplo is set to a LBFGSOperator of suitable size, but it may change to a LSR1Operator later on.

source
PartitionedStructures.Link.eplo_lsr1_from_epvMethod
eplo = eplo_lsr1_from_epv(epv::T) where {Y <: Number, T <: Elemental_pv{Y}}

Create an elemental limited-memory partitioned quasi-Newton operator PLSR1 eplo with the same partitioned structure than epv. Each element linear-operator of eplo is set to a LSR1Operator of suitable size.

source
PartitionedStructures.Link.epm_from_epvMethod
epm = epm_from_epv(epv::T) where {Y <: Number, T <: Elemental_pv{Y}}

Create an elemental partitioned quasi-Newton operator epm with the same partitioned structure than epv. Each element-matrix of epm is set with an identity matrix of suitable size.

source
PartitionedStructures.Link.epv_from_eploMethod
epv = epv_from_eplo(eplo::T) where T <: Part_mat{Y}

Create an elemental partitioned-vector epv with the same partitioned structure than eplo. Each element-vector of epv is set to a random vector of suitable size. Make a name difference with the method epv_from_epm().

source
PartitionedStructures.Link.epv_from_epmMethod
epv = epv_from_epm(epm::T) where T <: Part_mat{Y}

Create an elemental partitioned-vector epv with the same partitioned structure than epm. Each element-vector of epv is set to a random vector of suitable size.

source
PartitionedStructures.Link.mul_epm_epv!Method
mul_epm_epv!(epv_res::Elemental_pv{Y}, epm::T, epv::Elemental_pv{Y}) where {Y <: Number, T <: Part_mat{Y}}

Compute the elementwise product between the elemental partitioned-matrix epm and the elemental partitioned-vector epv. The result of each element-matrix element-vector product is stored in the elemental partitioned-vector epv_res.

source
PartitionedStructures.Link.mul_epm_epvMethod
epv_res = mul_epm_epv(epm::T, epv::Elemental_pv{Y}) where {Y <: Number, T <: Part_mat{Y}}

Compute the elementwise product between the elemental partitioned-matrix epm and the elemental partitioned-vector epv. The result is an elemental partitioned-vector epv_res storing the elementwise products between epm and epv.

source
PartitionedStructures.Link.mul_epm_vector!Method
mul_epm_vector!(res::Vector{Y}, epm::T, x::Vector{Y}) where {Y <: Number, T <: Part_mat{Y}}
mul_epm_vector!(res::Vector{Y}, epm::T, epv::Elemental_pv{Y}, x::Vector{Y}) where {Y <: Number, T <: Part_mat{Y}}

Compute the product between the elemental partitioned-matrix epm and the vector x. The method uses temporary the elemental partitioned-vector epv. The result is stored in res, a vector similar to x.

source
PartitionedStructures.Link.mul_epm_vectorMethod
result = mul_epm_vector(epm::T, x::Vector{Y}) where {Y <: Number, T <: Part_mat{Y}}
result = mul_epm_vector(epm::T, epv::Elemental_pv{Y}, x::Vector{Y}) where {Y <: Number, T <: Part_mat{Y}}

Compute the product between the elemental partitioned-matrix epm <: Part_mat and the vector x. The method uses temporary the elemental partitioned-vector epv. The method returns result, a vector similar to x.

source
PartitionedStructures.Link.string_counters_iterMethod
s = string_counters_iter(pm::T; name = :PQN) where {T <: Part_mat}

Produce s::String that summarizes the partitioned update applied onto pm at the last iterate. The method accumulates the informations gathered by each element-counter during the last iterate.

source
PartitionedStructures.Link.string_counters_totalMethod
s = string_counters_total(pm::T; name = :PQN) where {T <: Part_mat}

Produce s::String that summarizes the partitioned update applied onto pm since its allocations. The method accumulates the informations gathered by each element-counter since their allocations.

source
PartitionedStructures.PartitionedQuasiNewton.PBFGS_update!Method
PBFGS_update!(epm_B::Elemental_pm{T}, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where T
PBFGS_update!(epm_B::Elemental_pm{T}, epv_y::Elemental_pv{T}, epv_s::Elemental_pv{T}; verbose=true, kwargs...) where T

Perform the PBFGS update onto the elemental partitioned-matrix epm_B, given the step s (or the element-steps epv_s) and the difference of elemental partitioned-gradients epv_y.

source
PartitionedStructures.PartitionedQuasiNewton.PBFGS_updateMethod
copy_epm_B = PBFGS_update(epm_B::Elemental_pm{T}, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where T

Perform the PBFGS update onto a copy of the elemental partitioned-matrix epm_B, given the step s and the difference of elemental partitioned-gradients epv_y. Return the updated copy of epm_B.

source
PartitionedStructures.PartitionedQuasiNewton.PCS_update!Method
PCS_update!(epm_B::Elemental_pm{T}, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where T
PCS_update!(epm_B::Elemental_pm{T}, epv_y::Elemental_pv{T}, epv_s::Elemental_pv{T}; verbose=true, kwargs...) where T

Perform the PCS update onto the elemental partitioned-matrix epm_B, given the step s (or the element-steps epv_s) and the difference of elemental partitioned-gradients epv_y.

source
PartitionedStructures.PartitionedQuasiNewton.PCS_updateMethod
B = PCS_update(epm_B::Elemental_pm{T}, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where T

Perform the PCS update onto a copy of the elemental partitioned-matrix epm_B, given the step s and the difference of elemental partitioned-gradients epv_y. Each elemental element-matrix eem is update given its convex field, if convex==true the elemental element-matrix eem is update with BFGS otherwise it is SR1.

source
PartitionedStructures.PartitionedQuasiNewton.PSE_update!Method
PSE_update!(epm_B::Elemental_pm{T}, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where T
PSE_update!(epm_B::Elemental_pm{T}, epv_y::Elemental_pv{T}, epv_s::Elemental_pv{T}; verbose=true, kwargs...) where T

Perform the PSE update onto the elemental partitioned-matrix epm_B, given the step s (or the element-steps epv_s) and the difference of elemental partitioned-gradients epv_y.

source
PartitionedStructures.PartitionedQuasiNewton.PSE_updateMethod
copy_epm_B = PSE_update(epm_B::Elemental_pm{T}, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where T

Perform the PSE update onto a copy of the elemental partitioned-matrix epm_B, given the step s and the difference of elemental partitioned-gradients epv_y. Return the updated copy of epm_B.

source
PartitionedStructures.PartitionedQuasiNewton.PSR1_update!Method
PSR1_update!(epm_B::Elemental_pm{T}, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where T
PSR1_update!(epm_B::Elemental_pm{T}, epv_y::Elemental_pv{T}, epv_s::Elemental_pv{T}; verbose=true, kwargs...) where T

Performs the PSR1 update onto the elemental partitioned-matrix epm_B, given the step s (or the element-steps epv_s) and the difference of elemental partitioned-gradients epv_y.

source
PartitionedStructures.PartitionedQuasiNewton.PSR1_updateMethod
copy_epm_B = PSR1_update(epm_B::Elemental_pm{T}, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where T

Perform the PSR1 update onto a copy of the elemental partitioned-matrix epm_B, given the step s and the difference of elemental partitioned-gradients epv_y. Return the updated copy of epm_B.

source
PartitionedStructures.PartitionedLOQuasiNewton.PLBFGS_update!Method
PLBFGS_update!(eplo_B::Elemental_plo_bfgs{T}, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where T
PLBFGS_update!(eplo_B::Elemental_plo_bfgs{T}, epv_y::Elemental_pv{T}, epv_s::Elemental_pv{T}; verbose=true, reset=true, kwargs...) where T

Perform the PLBFGS update onto the partitioned limited-memory operator eplo_B, given the step s (or the element-steps epv_s) and the difference of elemental partitioned-gradients epv_y.

source
PartitionedStructures.PartitionedLOQuasiNewton.PLBFGS_updateMethod
copy_eplo_B = PLBFGS_update(eplo_B::Elemental_plo_bfgs{T}, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where T

Perform the PLBFGS update onto a copy of the partitioned limited-memory operator eplo_B, given the step s and the difference of elemental partitioned-gradients epv_y. Return the updated copy of eplo_B.

source
PartitionedStructures.PartitionedLOQuasiNewton.PLSE_update!Method
PLSE_update!(eplo_B::Y, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where {T, Y <: Part_LO_mat{T}}
PLSE_update!(eplo_B::Y, epv_y::Elemental_pv{T}, epv_s::Elemental_pv{T}; ω = 1e-6, verbose=true, reset=4, kwargs...) where {T, Y <: Part_LO_mat{T}}

Perform the partitionned update PLSE onto the partitioned limited-memory operator eplo_B, given the step s (or the element-steps epv_s) and the difference of elemental partitioned-gradients epv_y. Each element linear-operator from eplo_B is either a LBFGSOperator or LSR1Operator. The update tries to apply a LBFGS update to every Bᵢ, but if the curvature condition yᵢᵀUᵢs > 0 is not satisfied it replaces the LBFGSOperator by a LSR1Operator and applies a LSR1 update. If Bᵢ is initally a LSR1Opeartor, we replace it by a LBFGSOperator if the curvature condition yᵢᵀUᵢs > 0 holds and we update it, otherwise the LSR1Operator Bᵢ is update.

source
PartitionedStructures.PartitionedLOQuasiNewton.PLSE_updateMethod
copy_eplo_B = PLSE_update(eplo_B::Y, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where {T, Y <: Part_LO_mat{T}}

Perform the partitionned update PLSE onto a copy of the partitioned limited-memory operator eplo_B, given the step s and the difference of elemental partitioned-gradients epv_y. Each element linear-operator from eplo_B is either a LBFGSOperator or LSR1Operator. The update tries to apply a LBFGS update to every Bᵢ, but if the curvature condition yᵢᵀUᵢs > 0 is not satisfied it replaces the LBFGSOperator by a LSR1Operator and applies a LSR1 update. If Bᵢ is initally a LSR1Opeartor, we replace it by a LBFGSOperator if the curvature condition yᵢᵀUᵢs > 0 holds and we update it, otherwise the LSR1Operator Bᵢ is update. Return the updated copy of eplo_B.

source
PartitionedStructures.PartitionedLOQuasiNewton.PLSR1_update!Method
PLSR1_update!(eplo_B::Elemental_plo_sr1{T}, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where T
PLSR1_update!(eplo_B::Elemental_plo_sr1{T}, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where T

Perform the PLSR1 update onto the partitioned limited-memory operator eplo_B, given the step s (or the element-steps epv_s) and the difference of elemental partitioned-gradients epv_y.

source
PartitionedStructures.PartitionedLOQuasiNewton.PLSR1_updateMethod
copy_eplo_B = PLSR1_update(eplo_B::Elemental_plo_sr1{T}, epv_y::Elemental_pv{T}, s::Vector{T}; kwargs...) where T

Perform the PSR1 update onto a copy of the partitioned limited-memory operator eplo_B, given the step s and the difference of elemental partitioned-gradients epv_y. Return the updated copy of eplo_B.

source
PartitionedStructures.PartitionedLOQuasiNewton.Part_update!Method
Part_update!(eplo_B::Y, epv_y::Elemental_pv{T}, s::Vector{T}) where {T, Y <: Part_LO_mat{T}}
Part_update!(eplo_B::Y, epv_y::Elemental_pv{T}, epv_s::Elemental_pv{T}; kwargs...) where {T, Y <: Part_LO_mat{T}}

Perform a partitioned quasi-Newton update onto the partitioned limited-memory operator eplo_B, given the step s (or the element-steps epv_s) and the difference of elemental partitioned-gradients epv_y. Each elemental element linear-operator from eplo_B is either a LBFGSOperator or LSR1Operator. The update performs on each element the quasi-Newton update associated to the linear-operator.

source
PartitionedStructures.PartitionedLOQuasiNewton.Part_updateMethod
copy_eplo_B = Part_update(eplo_B::Y, epv_y::Elemental_pv{T}, s::Vector{T}) where {T, Y <: Part_LO_mat{T}}

Perform a quasi-Newton partitionned update onto a copy of the partitioned limited-memory operator eplo_B, given the step s and the difference of elemental partitioned-gradients epv_y. Each elemental element linear-operator from eplo_B is either a LBFGSOperator or LSR1Operator. The update performs on each element the quasi-Newton update associated to the linear-operator. Return the updated copy of eplo_B.

source
PartitionedStructures.ModElemental_elo_bfgs.Elemental_elo_bfgsType
Elemental_elo_bfgs{T} <: LOEltMat{T}

Represent an elemental element LBFGSOperator:

  • indices retains the indices of the elemental variables;
  • nie is the elemental size (=length(indices));
  • Bie a LBFGSOperator;
  • linear: if linear==true, then the element matrix contribution is null;
  • counter counts how many update the elemental limited-memory operator goes through from its allocation.
source
PartitionedStructures.ModElemental_plo_bfgs.Elemental_plo_bfgsType
Elemental_plo_bfgs{T} <: Part_LO_mat{T}

Represent an elemental partitioned quasi-Newton limited-memory operator PLBFGS. Each element is an elemental element LBFGSOperator. N is the number of elements. n is the size of the elemental partitioned limited-memory operator. eelo_set is the set of elemental element linear-operators. spm and L are sparse matrices either to form the sparse matrix gathering the elements or the Cholesky factor of spm. component_list summarizes for each variable i (∈ {1,..., n}) the list of elements (⊆ {1,...,N}) being parametrised by i. permutation is the current permutation of the elemental partitioned limited-memory operator ([1:n;] initially).

source
PartitionedStructures.ModElemental_plo_bfgs.PLBFGS_eploMethod
eplo = PLBFGS_eplo(; n::Int=9, T=Float64, nie::Int=5, overlapping::Int=1)

Return an elemental partitioned limited-memory operator PLBFGS of N (deduced from n and nie) elemental element linear-operators. Each element overlaps the coordinates of the next element by overlapping components.

source
PartitionedStructures.ModElemental_plo_bfgs.identity_eplo_LBFGSMethod
eplo = identity_eplo_LBFGS(element_variables::Vector{Vector{Int}}; N::Int=length(element_variables), n::Int=max_indices(element_variables), T=Float64, linear_vector::Vector{Bool} = zeros(Bool, N), mem=5)
eplo = identity_eplo_LBFGS(element_variables::Vector{Vector{Int}}, N::Int, n::Int; T=Float64, linear_vector::Vector{Bool} = zeros(Bool, N), mem=5)

Return an elemental partitioned limited-memory operator PLBFGS of N elemental element linear-operators. The positions are given by the vector of the element variables element_variables. linear_vector indicates which element linear-opeartor should not contribute to the partitioned linear-operator.

source
PartitionedStructures.ModElemental_elo_sr1.Elemental_elo_sr1Type
Elemental_elo_sr1{T} <: LOEltMat{T}

Represent an elemental element LSR1Operator;

  • indices retains the indices of the elemental variables;
  • nie is the elemental size (=length(indices));
  • Bie a LSR1Operator;
  • linear: if linear==true, then the element matrix contribution is null;
  • counter counts how many update the elemental limited-memory operator goes through from its allocation.
source
PartitionedStructures.ModElemental_plo_sr1.Elemental_plo_sr1Type
Elemental_plo_sr1{T} <: Part_LO_mat{T}

Represent an elemental partitioned quasi-Newton limited-memory operator PLSR1. Each element is an elemental element LSR1Operator. N is the number of elements. n is the size of the elemental partitioned limited-memory operator. eelo_set is the set of elemental element linear-operators. spm and L are sparse matrices either to form the sparse matrix gathering the elements or the Cholesky factor of spm. component_list summarizes for each variable i (∈ {1,..., n}) the list of elements (⊆ {1,...,N}) being parametrised by i. permutation is the current permutation of the elemental partitioned limited-memory operator ([1:n;] initially).

source
PartitionedStructures.ModElemental_plo_sr1.PLSR1_eploMethod
eplo = PLSR1_eplo(; n::Int=9, T=Float64, nie::Int=5, overlapping::Int=1)

Return an elemental partitionned limited-memory operator PLSR1 of N (deduced from n and nie) elemental element linear-operators. Each element overlaps the coordinates of the next element by overlapping components.

source
PartitionedStructures.ModElemental_plo_sr1.identity_eplo_LSR1Method
eplo = identity_eplo_LSR1(element_variables::Vector{Vector{Int}}; N::Int=length(element_variables), n::Int=max_indices(element_variables), T=Float64, linear_vector::Vector{Bool} = zeros(Bool, N); mem=5)
eplo = identity_eplo_LSR1(element_variables::Vector{Vector{Int}}, N::Int, n::Int; T=Float64, linear_vector::Vector{Bool} = zeros(Bool, N); mem=5)

Return an elemental partitionned limited-memory operator PLSR1 of N elemental element linear-operators. The positions are given by the vector of the element variables element_variables. linear_vector indicates which element linear-opeartor should not contribute to the partitioned linear-operator.

source
PartitionedStructures.ModElemental_plo.Elemental_ploType
Elemental_plo{T} <: Part_LO_mat{T}

Represent an elemental partitioned quasi-Newton limited-memory operator PLSE. Each element may either be a LBFGSOperator or a LSR1Operator. N is the number of elements. n is the size of the elemental partitioned limited-memory operator. eelo_set is the set of elemental element linear-operators. spm and L are sparse matrices either to form the sparse matrix gathering the elements or the Cholesky factor of spm. component_list summarizes for each variable i (∈ {1,..., n}) the list of elements (⊆ {1,...,N}) being parametrised by i. permutation is the current permutation of the elemental partitioned limited-memory operator ([1:n;] initially).

source
PartitionedStructures.ModElemental_plo.PLBFGSR1_eploMethod
eplo = PLBFGSR1_eplo(;n::Int=9, T=Float64, nie::Int=5, overlapping::Int=1, prob=0.5)

Create an elemental partitionned limited-memory operator PLSE of N (deduced from n and nie) elemental element linear-operators. Each element overlaps the coordinates of the next element by overlapping components. Each element is randomly (rand() > p) choose between an elemental element LBFGS operator or an elemental element LSR1 operator.

source
PartitionedStructures.ModElemental_plo.PLBFGSR1_eplo_randMethod
eplo = PLBFGSR1_eplo_rand(N::Int, n ::Int; T=Float64, nie::Int=5, prob=0.5)

Create an elemental partitionned limited-memory operator PLSE of N elemental element linear-operators. The size of each element is nie, whose positions are random in the range 1:n. Each element is randomly (rand() > p) choose between an elemental element LBFGS operator or an elemental element LSR1 operator.

source
PartitionedStructures.ModElemental_plo.identity_eplo_LOSEMethod
eplo = identity_eplo_LOSE(element_variables::Vector{Vector{Int}}; N::Int=length(element_variables), n::Int=max_indices(element_variables), T=Float64, linear_vector::Vector{Bool} = zeros(Bool, N), mem=5)
eplo = identity_eplo_LOSE(element_variables::Vector{Vector{Int}}, N::Int, n::Int; T=Float64, linear_vector::Vector{Bool} = zeros(Bool, N), mem=5)

Create an elemental partitionned limited-memory operator of N elemental element linear-operators initialized with LBFGS operators. The positions are given by the vector of the element variables element_variables. linear_vector indicates which element linear-opeartor should not contribute to the partitioned linear-operator.

source
PartitionedStructures.Instances.create_epv_eploMethod
(eplo,epv) = create_epv_eplo_sr1(;n=9,nie=5,overlapping=1, mul_v=100.)

Create an elemental partitioned limited-memory quasi-Newton operator eplo and elemental partitioned-vector epv. Both have the same partitioned structure defined by the size of the problem n::Int, the size of the element nie::Int and the overlapping between the consecutive elements overlapping::Int. Each elemental element-matrix is instantiated as a LBFGSOperator, but it may change to a LSR1Operator later on. The value of each elemental element-vector is rand(nie) .* mul_v::Real. Warning: You have to choose carefully the values n, nie and overlap, otherwise the method may fail. The default values are correct.

source
PartitionedStructures.Instances.create_epv_eplo_bfgsMethod
(eplo,epv) = create_epv_eplo_bfgs(;n=9,nie=5,overlapping=1, mul_v=100.)

Create an elemental partitioned limited-memory quasi-Newton operator PLBFGS eplo and an elemental partitioned-vector epv. Both have the same partitioned structure defined by the size of the problem n::Int, the size of the element nie::Int and the overlapping between the consecutive elements overlapping::Int. Each elemental element-matrix is a LBFGSOperator. The value of each elemental element-vector is rand(nie) .* mul_v::Real. Warning: You have to choose carefully the values n, nie and overlap, otherwise the method may fail. The default values are correct.

source
PartitionedStructures.Instances.create_epv_eplo_sr1Method
(eplo,epv) = create_epv_eplo_sr1(;n=9,nie=5,overlapping=1, mul_v=100.)

Create an elemental partitioned limited-memory quasi-Newton operator PLSR1 eplo and an elemental partitioned-vector epv. Both have the same partitioned structure defined by the size of the problem n::Int, the size of the element nie::Int and the overlapping between the consecutive elements overlapping::Int. Each elemental element-matrix is a LSR1Operator. The value of each elemental element-vector is rand(nie) .* mul_v::Real. Warning: You have to choose carefully the values n, nie and overlap, otherwise the method may fail. The default values are correct.

source
PartitionedStructures.Instances.create_epv_epmMethod
(epm,epv) = create_epv_epm(;n=9,nie=5,overlapping=1,mul_m=5., mul_v=100.)

Create an elemental partitioned-matrix epm and an elemental partitioned-vector epv. Both have the same partitioned structure defined by the size of the problem n::Int, the size of the element nie::Int and the overlapping between the consecutive elements overlapping::Int. Each elemental element-matrix is fill with ones, except the terms of the diagonal of value mul_v::Real. The value of each elemental element-vector is rand(nie) .* mul_v::Real. Warning: You have to choose carefully the values n, nie and overlap, otherwise the method may fail. The default values are correct.

source
PartitionedStructures.Instances.create_epv_epm_randMethod
(epm,epv) = create_epv_epm_rand(;n=9,nie=5,overlapping=1,range_mul_m=nie:2*nie, mul_v=100.)

Create an elemeental partitioned quasi-Newton operator epm and an elemental partitioned-vector epv. Both have the same partitioned structure defined by the size of the problem n::Int, the size of the element nie::Int and the overlapping between the consecutive elements overlapping::Int. Each elemental element-matrix is fill with ones, except the terms of the diagonal of value rand(1:range_mul_v). The value of each elemental element-vector is rand(nie) .* mul_v::Real. Warning: You have to choose carefully the values n, nie and overlap, otherwise the method may fail. The default values are correct.

source
PartitionedStructures.PartMatInterface.update!Method
update!(epm::T, epv_y::Elemental_pv{Y}, epv_s::Elemental_pv{Y}; name=:pse, kwargs...) where {Y <: Number, T <: Part_mat{Y}}
update!(epm::T, epv_y::Elemental_pv{Y}, s::Vector{Y}; kwargs...) where {Y <: Number, T <: Part_mat{Y}}

Update the elemental partitioned-matrix epm with a partitioned quasi-Newton update considering the difference of elemental partitioned-gradients epv_y and the step s (or elemental steps epv_s). The PSE update is run by default, you can apply a PBFGS or a PSR1 update with the optionnal argument name, respectively name=:pbfgs or name=:psr1.

source
PartitionedStructures.PartMatInterface.update!Method
update!(eplo::Elemental_plo_bfgs{Y}, epv_y::Elemental_pv{Y}, s::Vector{Y}; kwargs...) where Y <: Number
update!(eplo::Elemental_plo_bfgs{Y}, epv_y::Elemental_pv{Y}, epv_s::Elemental_pv{Y}; kwargs...) where Y <: Number

Updates the elemental partitioned limited-memory operator eplo with the partitioned quasi-Newton update PLBFGS considering the difference of elemental partitioned-gradients epv_y and the step s (or elemental steps epv_s).

source
PartitionedStructures.PartMatInterface.update!Method
update!(eplo::Elemental_plo_sr1{Y}, epv_y::Elemental_pv{Y}, s::Vector{Y}; kwargs...) where Y <: Number
update!(eplo::Elemental_plo_sr1{Y}, epv_y::Elemental_pv{Y}, epv_s::Elemental_pv{Y}; kwargs...) where Y <: Number

Updates the elemental partitioned limited-memory operator eplo with the partitioned quasi-Newton update PLSR1 considering the difference of elemental partitioned-gradients epv_y and the step s (or elemental steps epv_s).

source
PartitionedStructures.PartMatInterface.update!Method
update!(eplo::Elemental_plo{Y}, epv_y::Elemental_pv{Y}, s::Vector{Y}; kwargs...) where Y <: Number
update!(eplo::Elemental_plo{Y}, epv_y::Elemental_pv{Y}, epv_s::Elemental_pv{Y}; kwargs...) where Y <: Number

Updates the elemental partitioned limited-memory operator eplo with the partitioned quasi-Newton update PLSE considering the difference of elemental partitioned-gradients epv_y and the step s (or elemental steps epv_s).

source
PartitionedStructures.PartMatInterface.updateMethod
B = update(epm::T, epv_y::Elemental_pv{Y}, s::Vector{Y}; kwargs...) where {Y <: Number, T <: Part_mat{Y}}

Update the elemental partitioned operator epm <: Part_mat with a partitioned quasi-Newton update considering the difference of elemental partitioned-gradients epv_y and the step s. If epm is an elemental partitioned-matrix, the PSE update is run by default. You can apply a PBFGS or a PSR1 update with the optionnal argument name, respectively name=:pbfgs or name=:psr1. It returns a matrix accumulating every element-contribtion of the updated epm. Warning: this method should be use to test your algorithm, if you don't intend to form the matrix use update!(epm, epv_y, s).

source