Skip to content

Commit

Permalink
Merge pull request #49 from JuliaParallel/pr48
Browse files Browse the repository at this point in the history
update libelemental to v0.87.7 v2
  • Loading branch information
andreasnoack authored Jan 31, 2019
2 parents 94e1a5e + 44a65f6 commit 71db6e9
Show file tree
Hide file tree
Showing 4 changed files with 43 additions and 15 deletions.
5 changes: 3 additions & 2 deletions deps/build.jl
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
import Libdl, LibGit2, LinearAlgebra

# Use this version of Elemental
Elsha = "79987d38b04838acf6b6195be1967177521ee908"
# Use Elemental version 0.87.7
Elsha = "477e503a7a840cc1a75173552711b980505a0b06"

if Sys.iswindows()
error("Elemental only works on Unix Platforms")
Expand Down Expand Up @@ -49,6 +49,7 @@ cd(builddir) do
-D MATH_LIBS=$mathlib
-D EL_BLAS_SUFFIX=$blas_suffix
-D EL_LAPACK_SUFFIX=$blas_suffix
-D CMAKE_INSTALL_LIBDIR=$prefix/lib
-D CMAKE_INSTALL_RPATH=$prefix/lib
$srcdir`)
run(`make -j $build_procs`)
Expand Down
2 changes: 1 addition & 1 deletion src/io.jl
Original file line number Diff line number Diff line change
Expand Up @@ -4,7 +4,7 @@ for (elty, ext) in ((:Float32, :s),
(:Float64, :d),
(:ComplexF32, :c),
(:ComplexF64, :z))
for mattype in ("", "Dist")
for mattype in ("", "Dist", "Sparse", "DistSparse")
mat = Symbol(mattype, "Matrix")
@eval begin
function print(A::$mat{$elty}, title::String)
Expand Down
25 changes: 23 additions & 2 deletions src/optimization/solvers.jl
Original file line number Diff line number Diff line change
Expand Up @@ -33,6 +33,13 @@ struct MehrotraCtrl{T<:ElFloatType}
ruizMaxIter::ElInt
diagEquilTol::T
checkResiduals::ElBool

reg0Tmp::T
reg1Tmp::T
reg2Tmp::T
reg0Perm::T
reg1Perm::T
reg2Perm::T
end
function MehrotraCtrl(::Type{T};
primalInit::Bool = false,
Expand All @@ -55,7 +62,15 @@ function MehrotraCtrl(::Type{T};
ruizEquilTol = eps(T)^(-0.25),
ruizMaxIter = 3,
diagEquilTol = eps(T)^(-0.15),
checkResiduals = false) where {T<:ElFloatType}
checkResiduals = false,

reg0Tmp = eps(T)^0.25,
reg1Tmp = eps(T)^0.25,
reg2Tmp = eps(T)^0.25,

reg0Perm = eps(T)^0.35,
reg1Perm = eps(T)^0.35,
reg2Perm = eps(T)^0.35) where {T<:ElFloatType}

MehrotraCtrl{T}(ElBool(primalInit),
ElBool(dualInit),
Expand All @@ -77,7 +92,13 @@ function MehrotraCtrl(::Type{T};
T(ruizEquilTol),
ElInt(ruizMaxIter),
T(diagEquilTol),
ElBool(checkResiduals))
ElBool(checkResiduals),
T(reg0Tmp),
T(reg1Tmp),
T(reg2Tmp),
T(reg0Perm),
T(reg1Perm),
T(reg2Perm))
end

struct LPAffineCtrl{T<:ElFloatType}
Expand Down
26 changes: 16 additions & 10 deletions test/lav.jl
Original file line number Diff line number Diff line change
Expand Up @@ -22,6 +22,8 @@ function randDArray(m,n,sparsity=0.01)
return A
end

# Stack two 2D finite-difference matrices on top of each other
# and make the last column dense
function stackedFD2D(n0, n1)
height = 2*n0*n1
width = n0*n1
Expand Down Expand Up @@ -67,7 +69,7 @@ function stackedFD2D(n0, n1)
end

# The dense last column
El.queueUpdate(A, s, width, -div(10.0, height))
El.queueUpdate(A, s, width, floor(-10/height))
end
El.processQueues(A)
return A
Expand All @@ -88,16 +90,20 @@ bHeight = El.height(b)
bWidth = El.width(b)

if display
show(IO, A)
El.print(A, "A")
# El.print(b, "b")
end
ctrl = El.LPAffineCtrl(Float64,
mehrotraCtrl = El.MehrotraCtrl(Float64,
solveCtrl = El.RegSolveCtrl(Float64, progress = true),
print = true,
outerEquil = true,
time = true))

elapsedLAV = @elapsed x = El.lav(A, b)
ctrl = El.LPAffineCtrl(Float64,
mehrotraCtrl=El.MehrotraCtrl(Float64,
solveCtrl=El.RegSolveCtrl(Float64,
progress=true),
print=true,
outerEquil=true,
time=true)
)

# elapsedLAV = @elapsed x = El.lav(A, b);#Elemental.print(A, "Matrix A")
elapsedLAV = @elapsed x = El.lav(A, b, ctrl)

if El.MPI.commRank(El.MPI.CommWorld[]) == 0
Expand Down Expand Up @@ -128,7 +134,7 @@ end
rLS = copy(b)
mul!(rLS, A, xLS, -1.0, 1.0)
if display
El.Display( rLS, "A x_{LS} - b" )
El.print( rLS, "A x_{LS} - b" )
end

rLSTwoNorm = El.nrm2(rLS)
Expand Down

0 comments on commit 71db6e9

Please sign in to comment.