1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501
|
function spqr_make (opt1)
%SPQR_MAKE compiles the SuiteSparseQR mexFunctions
%
% Example:
% spqr_make
%
% SuiteSparseQR relies on CHOLMOD, AMD, and COLAMD, and optionally CCOLAMD,
% CAMD, and METIS. Next, type
%
% spqr_make
%
% in the MATLAB command window. If METIS is not present in ../../metis-5.1.0,
% then it is not used.
%
% To compile using Intel's Threading Building Blocks (TBB) use:
%
% spqr_make ('tbb')
%
% TBB parallelism is not the default, since it conflicts with the multithreaded
% BLAS (the Intel MKL are OpenMP based, for example). This may change in
% future versions.
%
% You must type the spqr_make command while in the SuiteSparseQR/MATLAB
% directory.
%
% See also spqr, spqr_solve, spqr_qmult, qr, mldivide
% Copyright 2008, Timothy A. Davis, http://www.suitesparse.com
details = 0 ; % 1 if details of each command are to be printed, 0 if not
v = version ;
try
% ispc does not appear in MATLAB 5.3
pc = ispc ;
mac = ismac ;
catch %#ok
% if ispc fails, assume we are on a Windows PC if it's not unix
pc = ~isunix ;
mac = 0 ;
end
flags = '' ;
is64 = (~isempty (strfind (computer, '64'))) ;
if (is64)
% 64-bit MATLAB
flags = '-largeArrayDims' ;
end
% MATLAB 8.3.0 now has a -silent option to keep 'mex' from burbling too much
if (~verLessThan ('matlab', '8.3.0'))
flags = ['-silent ' flags] ;
end
include = '-DNMATRIXOPS -DNMODIFY -I. -I../../AMD/Include -I../../COLAMD/Include -I../../CHOLMOD/Include -I../Include -I../../SuiteSparse_config' ;
% Determine if METIS is available
metis_path = '../../metis-5.1.0' ;
have_metis = exist (metis_path, 'dir') ;
% Determine if TBB is to be used
if (nargin < 1)
tbb = 0 ;
elseif (nargin < 2)
tbb = strcmp (opt1, 'tbb') ;
end
% fix the METIS 4.0.1 rename.h file
if (have_metis)
fprintf ('Compiling SuiteSparseQR with METIS for MATLAB Version %s\n', v) ;
include = [include ' -I' metis_path '/include'] ;
include = [include ' -I' metis_path '/GKlib'] ;
include = [include ' -I' metis_path '/libmetis'] ;
include = [include ' -I../../CCOLAMD/Include -I../../CAMD/Include' ] ;
else
fprintf ('Compiling SuiteSparseQR without METIS on MATLAB Version %s\n', v);
include = ['-DNPARTITION ' include ] ;
end
%-------------------------------------------------------------------------------
% BLAS option
%-------------------------------------------------------------------------------
% This is exceedingly ugly. The MATLAB mex command needs to be told where to
% find the LAPACK and BLAS libraries, which is a real portability nightmare.
% The correct option is highly variable and depends on the MATLAB version.
if (pc)
if (verLessThan ('matlab', '6.5'))
% MATLAB 6.1 and earlier: use the version supplied in CHOLMOD
lib = '../../CHOLMOD/MATLAB/lcc_lib/libmwlapack.lib' ;
elseif (verLessThan ('matlab', '7.5'))
% use the built-in LAPACK lib (which includes the BLAS)
lib = 'libmwlapack.lib' ;
elseif (verLessThan ('matlab', '9.5'))
lib = 'libmwlapack.lib libmwblas.lib' ;
else
lib = '-lmwlapack -lmwblas' ;
end
else
if (verLessThan ('matlab', '7.5'))
% MATLAB 7.5 and earlier, use the LAPACK lib (including the BLAS)
lib = '-lmwlapack' ;
else
% MATLAB 7.6 requires the -lmwblas option; earlier versions do not
lib = '-lmwlapack -lmwblas' ;
end
end
if (is64 && ~verLessThan ('matlab', '7.8'))
% versions 7.8 and later on 64-bit platforms use a 64-bit BLAS
fprintf ('with 64-bit BLAS\n') ;
flags = [flags ' -DBLAS64'] ;
end
%-------------------------------------------------------------------------------
% GPU option
%-------------------------------------------------------------------------------
% GPU not yet supported for the spqr MATLAB mexFunction
% flags = [flags ' -DGPU_BLAS'] ;
%-------------------------------------------------------------------------------
% TBB option
%-------------------------------------------------------------------------------
% You should install TBB properly so that mex can find the library files and
% include files, but you can also modify the tbb_lib_path and tbb_include_path
% strings below to if you need to specify the path to your own installation of
% TBB.
% vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv
% >>>>>>>>>>>>>>>>>>>>> EDIT THE tbb_path BELOW AS NEEDED <<<<<<<<<<<<<<<<<<<<<<
if (pc)
% For Windows, with TBB installed in C:\TBB. Edit this line as needed:
tbb_path = 'C:\TBB\tbb21_009oss' ;
else
% For Linux, edit this line as needed (not needed if already in /usr/lib):
tbb_path = '' ;
end
% ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
% You should not have to edit the lines below.
if (pc)
if (is64)
tbb_lib_path = [tbb_path '\ia32\vc9\lib\'] ;
else
tbb_lib_path = [tbb_path '\em64t\vc9\lib\'] ;
end
tbb_include_path = [tbb_path '\include\'] ;
else
% For Linux, with TBB might be already installed in /usr/lib
if (exist ('/usr/lib/libtbb.so', 'file'))
% do not edit these lines
tbb_path = '' ;
tbb_lib_path = '' ;
tbb_include_path = '' ;
else
if (is64)
tbb_lib_path = '/em64t/cc4.1.0_libc2.4_kernel2.6.16.21/lib' ;
else
tbb_lib_path = '/ia32/cc4.1.0_libc2.4_kernel2.6.16.21/lib' ;
end
tbb_lib_path = [tbb_path tbb_lib_path] ;
tbb_include_path = [tbb_path '/include/'] ;
end
end
if (tbb)
fprintf ('Compiling with Intel TBB parallelism\n') ;
lib = [lib ' -L' tbb_lib_path ' -ltbb'] ;
include = [include ' -I' tbb_include_path ' -DHAVE_TBB' ] ;
end
if (~(pc || mac))
% for POSIX timing routine
lib = [lib ' -lrt'] ;
end
%-------------------------------------------------------------------------------
% ready to compile ...
%-------------------------------------------------------------------------------
config_src = { '../../SuiteSparse_config/SuiteSparse_config' } ;
amd_c_src = { ...
'../../AMD/Source/amd_1', ...
'../../AMD/Source/amd_2', ...
'../../AMD/Source/amd_aat', ...
'../../AMD/Source/amd_control', ...
'../../AMD/Source/amd_defaults', ...
'../../AMD/Source/amd_dump', ...
'../../AMD/Source/amd_global', ...
'../../AMD/Source/amd_info', ...
'../../AMD/Source/amd_order', ...
'../../AMD/Source/amd_postorder', ...
'../../AMD/Source/amd_post_tree', ...
'../../AMD/Source/amd_preprocess', ...
'../../AMD/Source/amd_valid' } ;
colamd_c_src = {
'../../COLAMD/Source/colamd' } ;
% CAMD and CCOLAMD are not needed if we don't have METIS
camd_c_src = { ...
'../../CAMD/Source/camd_1', ...
'../../CAMD/Source/camd_2', ...
'../../CAMD/Source/camd_aat', ...
'../../CAMD/Source/camd_control', ...
'../../CAMD/Source/camd_defaults', ...
'../../CAMD/Source/camd_dump', ...
'../../CAMD/Source/camd_global', ...
'../../CAMD/Source/camd_info', ...
'../../CAMD/Source/camd_order', ...
'../../CAMD/Source/camd_postorder', ...
'../../CAMD/Source/camd_preprocess', ...
'../../CAMD/Source/camd_valid' } ;
ccolamd_c_src = {
'../../CCOLAMD/Source/ccolamd' } ;
if (have_metis)
metis_c_src = {
'GKlib/b64', ...
'GKlib/blas', ...
'GKlib/csr', ...
'GKlib/error', ...
'GKlib/evaluate', ...
'GKlib/fkvkselect', ...
'GKlib/fs', ...
'GKlib/getopt', ...
'GKlib/gkregex', ...
'GKlib/graph', ...
'GKlib/htable', ...
'GKlib/io', ...
'GKlib/itemsets', ...
'GKlib/mcore', ...
'GKlib/memory', ...
'GKlib/omp', ...
'GKlib/pdb', ...
'GKlib/pqueue', ...
'GKlib/random', ...
'GKlib/rw', ...
'GKlib/seq', ...
'GKlib/sort', ...
'GKlib/string', ...
'GKlib/timers', ...
'GKlib/tokenizer', ...
'GKlib/util', ...
'libmetis/auxapi', ...
'libmetis/balance', ...
'libmetis/bucketsort', ...
'libmetis/checkgraph', ...
'libmetis/coarsen', ...
'libmetis/compress', ...
'libmetis/contig', ...
'libmetis/debug', ...
'libmetis/fm', ...
'libmetis/fortran', ...
'libmetis/frename', ...
'libmetis/gklib', ...
'libmetis/graph', ...
'libmetis/initpart', ...
'libmetis/kmetis', ...
'libmetis/kwayfm', ...
'libmetis/kwayrefine', ...
'libmetis/mcutil', ...
'libmetis/mesh', ...
'libmetis/meshpart', ...
'libmetis/minconn', ...
'libmetis/mincover', ...
'libmetis/mmd', ...
'libmetis/ometis', ...
'libmetis/options', ...
'libmetis/parmetis', ...
'libmetis/pmetis', ...
'libmetis/refine', ...
'libmetis/separator', ...
'libmetis/sfm', ...
'libmetis/srefine', ...
'libmetis/stat', ...
'libmetis/timing', ...
'libmetis/util', ...
'libmetis/wspace', ...
} ;
for i = 1:length (metis_c_src)
metis_c_src {i} = [metis_path '/' metis_c_src{i}] ;
end
end
cholmod_c_src = {
'../../CHOLMOD/Core/cholmod_aat', ...
'../../CHOLMOD/Core/cholmod_add', ...
'../../CHOLMOD/Core/cholmod_band', ...
'../../CHOLMOD/Core/cholmod_change_factor', ...
'../../CHOLMOD/Core/cholmod_common', ...
'../../CHOLMOD/Core/cholmod_complex', ...
'../../CHOLMOD/Core/cholmod_copy', ...
'../../CHOLMOD/Core/cholmod_dense', ...
'../../CHOLMOD/Core/cholmod_error', ...
'../../CHOLMOD/Core/cholmod_factor', ...
'../../CHOLMOD/Core/cholmod_memory', ...
'../../CHOLMOD/Core/cholmod_sparse', ...
'../../CHOLMOD/Core/cholmod_transpose', ...
'../../CHOLMOD/Core/cholmod_triplet', ...
'../../CHOLMOD/Check/cholmod_check', ...
'../../CHOLMOD/Check/cholmod_read', ...
'../../CHOLMOD/Check/cholmod_write', ...
'../../CHOLMOD/Cholesky/cholmod_amd', ...
'../../CHOLMOD/Cholesky/cholmod_analyze', ...
'../../CHOLMOD/Cholesky/cholmod_colamd', ...
'../../CHOLMOD/Cholesky/cholmod_etree', ...
'../../CHOLMOD/Cholesky/cholmod_factorize', ...
'../../CHOLMOD/Cholesky/cholmod_postorder', ...
'../../CHOLMOD/Cholesky/cholmod_rcond', ...
'../../CHOLMOD/Cholesky/cholmod_resymbol', ...
'../../CHOLMOD/Cholesky/cholmod_rowcolcounts', ...
'../../CHOLMOD/Cholesky/cholmod_rowfac', ...
'../../CHOLMOD/Cholesky/cholmod_solve', ...
'../../CHOLMOD/Cholesky/cholmod_spsolve', ...
'../../CHOLMOD/Supernodal/cholmod_super_numeric', ...
'../../CHOLMOD/Supernodal/cholmod_super_solve', ...
'../../CHOLMOD/Supernodal/cholmod_super_symbolic' } ;
cholmod_c_partition_src = {
'../../CHOLMOD/Partition/cholmod_ccolamd', ...
'../../CHOLMOD/Partition/cholmod_csymamd', ...
'../../CHOLMOD/Partition/cholmod_camd', ...
'../../CHOLMOD/Partition/cholmod_metis', ...
'../../CHOLMOD/Partition/cholmod_nesdis' } ;
% SuiteSparseQR does not need the MatrixOps or Modify modules of CHOLMOD
% cholmod_unused = {
% '../../CHOLMOD/MatrixOps/cholmod_drop', ...
% '../../CHOLMOD/MatrixOps/cholmod_horzcat', ...
% '../../CHOLMOD/MatrixOps/cholmod_norm', ...
% '../../CHOLMOD/MatrixOps/cholmod_scale', ...
% '../../CHOLMOD/MatrixOps/cholmod_sdmult', ...
% '../../CHOLMOD/MatrixOps/cholmod_ssmult', ...
% '../../CHOLMOD/MatrixOps/cholmod_submatrix', ...
% '../../CHOLMOD/MatrixOps/cholmod_vertcat', ...
% '../../CHOLMOD/MatrixOps/cholmod_symmetry', ...
% '../../CHOLMOD/Modify/cholmod_rowadd', ...
% '../../CHOLMOD/Modify/cholmod_rowdel', ...
% '../../CHOLMOD/Modify/cholmod_updown' } ;
% SuiteSparseQR source code, and mex support file
spqr_cpp_src = {
'../Source/spqr_parallel', ...
'../Source/spqr_1colamd', ...
'../Source/spqr_1factor', ...
'../Source/spqr_1fixed', ...
'../Source/spqr_analyze', ...
'../Source/spqr_append', ...
'../Source/spqr_assemble', ...
'../Source/spqr_cpack', ...
'../Source/spqr_csize', ...
'../Source/spqr_cumsum', ...
'../Source/spqr_debug', ...
'../Source/spqr_factorize', ...
'../Source/spqr_fcsize', ...
'../Source/spqr_freefac', ...
'../Source/spqr_freenum', ...
'../Source/spqr_freesym', ...
'../Source/spqr_front', ...
'../Source/spqr_fsize', ...
'../Source/spqr_happly', ...
'../Source/spqr_happly_work', ...
'../Source/spqr_hpinv', ...
'../Source/spqr_kernel', ...
'../Source/spqr_larftb', ...
'../Source/spqr_panel', ...
'../Source/spqr_rconvert', ...
'../Source/spqr_rcount', ...
'../Source/spqr_rhpack', ...
'../Source/spqr_rmap', ...
'../Source/spqr_rsolve', ...
'../Source/spqr_shift', ...
'../Source/spqr_stranspose1', ...
'../Source/spqr_stranspose2', ...
'../Source/spqr_trapezoidal', ...
'../Source/spqr_type', ...
'../Source/spqr_tol', ...
'../Source/spqr_maxcolnorm', ...
'../Source/SuiteSparseQR_qmult', ...
'../Source/SuiteSparseQR', ...
'../Source/SuiteSparseQR_expert', ...
'../MATLAB/spqr_mx' } ;
% SuiteSparse C source code, for MATLAB error handling
spqr_c_mx_src = { '../MATLAB/spqr_mx_error' } ;
% SuiteSparseQR mexFunctions
spqr_mex_cpp_src = { 'spqr', 'spqr_qmult', 'spqr_solve', 'spqr_singletons' } ;
if (pc)
% Windows does not have drand48 and srand48, required by METIS. Use
% drand48 and srand48 in CHOLMOD/MATLAB/Windows/rand48.c instead.
% Also provide Windows with an empty <strings.h> include file.
obj_extension = '.obj' ;
cholmod_c_src = [cholmod_c_src {'../../CHOLMOD/MATLAB/Windows/rand48'}] ;
include = [include ' -I../../CHOLMOD/MATLAB/Windows'] ;
else
obj_extension = '.o' ;
end
% compile each library source file
obj = '' ;
c_source = [config_src amd_c_src colamd_c_src cholmod_c_src spqr_c_mx_src ] ;
if (have_metis)
c_source = [c_source cholmod_c_partition_src ccolamd_c_src ] ;
c_source = [c_source camd_c_src metis_c_src] ;
end
cpp_source = spqr_cpp_src ;
kk = 0 ;
for f = cpp_source
ff = f {1} ;
slash = strfind (ff, '/') ;
if (isempty (slash))
slash = 1 ;
else
slash = slash (end) + 1 ;
end
o = ff (slash:end) ;
obj = [obj ' ' o obj_extension] ; %#ok
s = sprintf ('mex %s -O %s -c %s.cpp', flags, include, ff) ;
kk = do_cmd (s, kk, details) ;
end
for f = c_source
ff = f {1} ;
if (isequal (ff, [metis_path '/GKlib/util']))
% special case, since a file with the same name also exists in libmetis
copyfile ([ff '.c'], 'GKlib_util.c', 'f') ;
ff = 'GKlib_util' ;
o = 'GKlib_util' ;
elseif (isequal (ff, [metis_path '/GKlib/graph']))
% special case, since a file with the same name also exist in libmetis
copyfile ([ff '.c'], 'GKlib_graph.c', 'f') ;
ff = 'GKlib_graph' ;
o = 'GKlib_graph' ;
else
slash = strfind (ff, '/') ;
if (isempty (slash))
slash = 1 ;
else
slash = slash (end) + 1 ;
end
o = ff (slash:end) ;
end
% fprintf ('%s\n', o) ;
o = [o obj_extension] ;
obj = [obj ' ' o] ; %#ok
s = sprintf ('mex %s -DDLONG -O %s -c %s.c', flags, include, ff) ;
kk = do_cmd (s, kk, details) ;
end
% compile each mexFunction
for f = spqr_mex_cpp_src
s = sprintf ('mex %s -O %s %s.cpp', flags, include, f{1}) ;
s = [s obj ' ' lib] ; %#ok
kk = do_cmd (s, kk, details) ;
end
% clean up
s = ['delete ' obj] ;
status = warning ('off', 'MATLAB:DELETE:FileNotFound') ;
delete rename.h
warning (status) ;
do_cmd (s, kk, details) ;
fprintf ('\nSuiteSparseQR successfully compiled\n') ;
% remove the renamed METIS files, if they exist
if (exist ('GKlib_util.c', 'file'))
delete ('GKlib_util.c') ;
end
if (exist ('GKlib_graph.c', 'file'))
delete ('GKlib_graph.c') ;
end
%-------------------------------------------------------------------------------
function kk = do_cmd (s, kk, details)
%DO_CMD evaluate a command, and either print it or print a "."
if (details)
fprintf ('%s\n', s) ;
else
if (mod (kk, 60) == 0)
fprintf ('\n') ;
end
kk = kk + 1 ;
fprintf ('.') ;
end
eval (s) ;
|