00001 /* ========================================================================== */ 00002 /* === Include/cholmod_core.h =============================================== */ 00003 /* ========================================================================== */ 00004 00005 /* ----------------------------------------------------------------------------- 00006 * CHOLMOD/Include/cholmod_core.h. 00007 * Copyright (C) 2005-2006, Univ. of Florida. Author: Timothy A. Davis 00008 * CHOLMOD/Include/cholmod_core.h is licensed under Version 2.1 of the GNU 00009 * Lesser General Public License. See lesser.txt for a text of the license. 00010 * CHOLMOD is also available under other licenses; contact authors for details. 00011 * http://www.cise.ufl.edu/research/sparse 00012 * -------------------------------------------------------------------------- */ 00013 00014 /* CHOLMOD Core module: basic CHOLMOD objects and routines. 00015 * Required by all CHOLMOD modules. Requires no other module or package. 00016 * 00017 * The CHOLMOD modules are: 00018 * 00019 * Core basic data structures and definitions 00020 * Check check/print the 5 CHOLMOD objects, & 3 types of integer vectors 00021 * Cholesky sparse Cholesky factorization 00022 * Modify sparse Cholesky update/downdate/row-add/row-delete 00023 * MatrixOps sparse matrix functions (add, multiply, norm, ...) 00024 * Supernodal supernodal sparse Cholesky factorization 00025 * Partition graph-partitioning based orderings 00026 * 00027 * The CHOLMOD objects: 00028 * -------------------- 00029 * 00030 * cholmod_common parameters, statistics, and workspace 00031 * cholmod_sparse a sparse matrix in compressed column form 00032 * cholmod_factor an LL' or LDL' factorization 00033 * cholmod_dense a dense matrix 00034 * cholmod_triplet a sparse matrix in "triplet" form 00035 * 00036 * The Core module described here defines the CHOLMOD data structures, and 00037 * basic operations on them. To create and solve a sparse linear system Ax=b, 00038 * the user must create A and b, populate them with values, and then pass them 00039 * to the routines in the CHOLMOD Cholesky module. There are two primary 00040 * methods for creating A: (1) allocate space for a column-oriented sparse 00041 * matrix and fill it with pattern and values, or (2) create a triplet form 00042 * matrix and convert it to a sparse matrix. The latter option is simpler. 00043 * 00044 * The matrices b and x are typically dense matrices, but can also be sparse. 00045 * You can allocate and free them as dense matrices with the 00046 * cholmod_allocate_dense and cholmod_free_dense routines. 00047 * 00048 * The cholmod_factor object contains the symbolic and numeric LL' or LDL' 00049 * factorization of sparse symmetric matrix. The matrix must be positive 00050 * definite for an LL' factorization. It need only be symmetric and have well- 00051 * conditioned leading submatrices for it to have an LDL' factorization 00052 * (CHOLMOD does not pivot for numerical stability). It is typically created 00053 * with the cholmod_factorize routine in the Cholesky module, but can also 00054 * be initialized to L=D=I in the Core module and then modified by the Modify 00055 * module. It must be freed with cholmod_free_factor, defined below. 00056 * 00057 * The Core routines for each object are described below. Each list is split 00058 * into two parts: the primary routines and secondary routines. 00059 * 00060 * ============================================================================ 00061 * === cholmod_common ========================================================= 00062 * ============================================================================ 00063 * 00064 * The Common object contains control parameters, statistics, and 00065 * You must call cholmod_start before calling any other CHOLMOD routine, and 00066 * must call cholmod_finish as your last call to CHOLMOD, with two exceptions: 00067 * you may call cholmod_print_common and cholmod_check_common in the Check 00068 * module after calling cholmod_finish. 00069 * 00070 * cholmod_start first call to CHOLMOD 00071 * cholmod_finish last call to CHOLMOD 00072 * ----------------------------- 00073 * cholmod_defaults restore default parameters 00074 * cholmod_maxrank maximum rank for update/downdate 00075 * cholmod_allocate_work allocate workspace in Common 00076 * cholmod_free_work free workspace in Common 00077 * cholmod_clear_flag clear Flag workspace in Common 00078 * cholmod_error called when CHOLMOD encounters an error 00079 * cholmod_dbound for internal use in CHOLMOD only 00080 * cholmod_hypot compute sqrt (x*x + y*y) accurately 00081 * cholmod_divcomplex complex division, c = a/b 00082 * 00083 * ============================================================================ 00084 * === cholmod_sparse ========================================================= 00085 * ============================================================================ 00086 * 00087 * A sparse matrix is held in compressed column form. In the basic type 00088 * ("packed", which corresponds to a MATLAB sparse matrix), an n-by-n matrix 00089 * with nz entries is held in three arrays: p of size n+1, i of size nz, and x 00090 * of size nz. Row indices of column j are held in i [p [j] ... p [j+1]-1] and 00091 * in the same locations in x. There may be no duplicate entries in a column. 00092 * Row indices in each column may be sorted or unsorted (CHOLMOD keeps track). 00093 * A->stype determines the storage mode: 0 if both upper/lower parts are stored, 00094 * -1 if A is symmetric and just tril(A) is stored, +1 if symmetric and triu(A) 00095 * is stored. 00096 * 00097 * cholmod_allocate_sparse allocate a sparse matrix 00098 * cholmod_free_sparse free a sparse matrix 00099 * ----------------------------- 00100 * cholmod_reallocate_sparse change the size (# entries) of sparse matrix 00101 * cholmod_nnz number of nonzeros in a sparse matrix 00102 * cholmod_speye sparse identity matrix 00103 * cholmod_spzeros sparse zero matrix 00104 * cholmod_transpose transpose a sparse matrix 00105 * cholmod_ptranspose transpose/permute a sparse matrix 00106 * cholmod_transpose_unsym transpose/permute an unsymmetric sparse matrix 00107 * cholmod_transpose_sym transpose/permute a symmetric sparse matrix 00108 * cholmod_sort sort row indices in each column of sparse matrix 00109 * cholmod_band C = tril (triu (A,k1), k2) 00110 * cholmod_band_inplace A = tril (triu (A,k1), k2) 00111 * cholmod_aat C = A*A' 00112 * cholmod_copy_sparse C = A, create an exact copy of a sparse matrix 00113 * cholmod_copy C = A, with possible change of stype 00114 * cholmod_add C = alpha*A + beta*B 00115 * cholmod_sparse_xtype change the xtype of a sparse matrix 00116 * 00117 * ============================================================================ 00118 * === cholmod_factor ========================================================= 00119 * ============================================================================ 00120 * 00121 * The data structure for an LL' or LDL' factorization is too complex to 00122 * describe in one sentence. This object can hold the symbolic analysis alone, 00123 * or in combination with a "simplicial" (similar to a sparse matrix) or 00124 * "supernodal" form of the numerical factorization. Only the routine to free 00125 * a factor is primary, since a factor object is created by the factorization 00126 * routine (cholmod_factorize). It must be freed with cholmod_free_factor. 00127 * 00128 * cholmod_free_factor free a factor 00129 * ----------------------------- 00130 * cholmod_allocate_factor allocate a factor (LL' or LDL') 00131 * cholmod_reallocate_factor change the # entries in a factor 00132 * cholmod_change_factor change the type of factor (e.g., LDL' to LL') 00133 * cholmod_pack_factor pack the columns of a factor 00134 * cholmod_reallocate_column resize a single column of a factor 00135 * cholmod_factor_to_sparse create a sparse matrix copy of a factor 00136 * cholmod_copy_factor create a copy of a factor 00137 * cholmod_factor_xtype change the xtype of a factor 00138 * 00139 * Note that there is no cholmod_sparse_to_factor routine to create a factor 00140 * as a copy of a sparse matrix. It could be done, after a fashion, but a 00141 * lower triangular sparse matrix would not necessarily have a chordal graph, 00142 * which would break the many CHOLMOD routines that rely on this property. 00143 * 00144 * ============================================================================ 00145 * === cholmod_dense ========================================================== 00146 * ============================================================================ 00147 * 00148 * The solve routines and some of the MatrixOps and Modify routines use dense 00149 * matrices as inputs. These are held in column-major order. With a leading 00150 * dimension of d, the entry in row i and column j is held in x [i+j*d]. 00151 * 00152 * cholmod_allocate_dense allocate a dense matrix 00153 * cholmod_free_dense free a dense matrix 00154 * ----------------------------- 00155 * cholmod_zeros allocate a dense matrix of all zeros 00156 * cholmod_ones allocate a dense matrix of all ones 00157 * cholmod_eye allocate a dense identity matrix 00158 * cholmod_sparse_to_dense create a dense matrix copy of a sparse matrix 00159 * cholmod_dense_to_sparse create a sparse matrix copy of a dense matrix 00160 * cholmod_copy_dense create a copy of a dense matrix 00161 * cholmod_copy_dense2 copy a dense matrix (pre-allocated) 00162 * cholmod_dense_xtype change the xtype of a dense matrix 00163 * 00164 * ============================================================================ 00165 * === cholmod_triplet ======================================================== 00166 * ============================================================================ 00167 * 00168 * A sparse matrix held in triplet form is the simplest one for a user to 00169 * create. It consists of a list of nz entries in arbitrary order, held in 00170 * three arrays: i, j, and x, each of length nk. The kth entry is in row i[k], 00171 * column j[k], with value x[k]. There may be duplicate values; if A(i,j) 00172 * appears more than once, its value is the sum of the entries with those row 00173 * and column indices. 00174 * 00175 * cholmod_allocate_triplet allocate a triplet matrix 00176 * cholmod_triplet_to_sparse create a sparse matrix copy of a triplet matrix 00177 * cholmod_free_triplet free a triplet matrix 00178 * ----------------------------- 00179 * cholmod_reallocate_triplet change the # of entries in a triplet matrix 00180 * cholmod_sparse_to_triplet create a triplet matrix copy of a sparse matrix 00181 * cholmod_copy_triplet create a copy of a triplet matrix 00182 * cholmod_triplet_xtype change the xtype of a triplet matrix 00183 * 00184 * ============================================================================ 00185 * === memory management ====================================================== 00186 * ============================================================================ 00187 * 00188 * cholmod_malloc malloc wrapper 00189 * cholmod_calloc calloc wrapper 00190 * cholmod_free free wrapper 00191 * cholmod_realloc realloc wrapper 00192 * cholmod_realloc_multiple realloc wrapper for multiple objects 00193 * 00194 * ============================================================================ 00195 * === Core CHOLMOD prototypes ================================================ 00196 * ============================================================================ 00197 * 00198 * All CHOLMOD routines (in all modules) use the following protocol for return 00199 * values, with one exception: 00200 * 00201 * int TRUE (1) if successful, or FALSE (0) otherwise. 00202 * (exception: cholmod_divcomplex) 00203 * UF_long a value >= 0 if successful, or -1 otherwise. 00204 * double a value >= 0 if successful, or -1 otherwise. 00205 * size_t a value > 0 if successful, or 0 otherwise. 00206 * void * a non-NULL pointer to newly allocated memory if 00207 * successful, or NULL otherwise. 00208 * cholmod_sparse * a non-NULL pointer to a newly allocated matrix 00209 * if successful, or NULL otherwise. 00210 * cholmod_factor * a non-NULL pointer to a newly allocated factor 00211 * if successful, or NULL otherwise. 00212 * cholmod_triplet * a non-NULL pointer to a newly allocated triplet 00213 * matrix if successful, or NULL otherwise. 00214 * cholmod_dense * a non-NULL pointer to a newly allocated triplet 00215 * matrix if successful, or NULL otherwise. 00216 * 00217 * The last parameter to all routines is always a pointer to the CHOLMOD 00218 * Common object. 00219 * 00220 * TRUE and FALSE are not defined here, since they may conflict with the user 00221 * program. A routine that described here returning TRUE or FALSE returns 1 00222 * or 0, respectively. Any TRUE/FALSE parameter is true if nonzero, false if 00223 * zero. 00224 */ 00225 00226 #ifndef CHOLMOD_CORE_H 00227 #define CHOLMOD_CORE_H 00228 00229 /* ========================================================================== */ 00230 /* === CHOLMOD version ====================================================== */ 00231 /* ========================================================================== */ 00232 00233 /* All versions of CHOLMOD will include the following definitions. 00234 * As an example, to test if the version you are using is 1.3 or later: 00235 * 00236 * if (CHOLMOD_VERSION >= CHOLMOD_VER_CODE (1,3)) ... 00237 * 00238 * This also works during compile-time: 00239 * 00240 * #if CHOLMOD_VERSION >= CHOLMOD_VER_CODE (1,3) 00241 * printf ("This is version 1.3 or later\n") ; 00242 * #else 00243 * printf ("This is version is earlier than 1.3\n") ; 00244 * #endif 00245 */ 00246 00247 #define CHOLMOD_DATE "Sept 30, 2008" 00248 #define CHOLMOD_VER_CODE(main,sub) ((main) * 1000 + (sub)) 00249 #define CHOLMOD_MAIN_VERSION 1 00250 #define CHOLMOD_SUB_VERSION 7 00251 #define CHOLMOD_SUBSUB_VERSION 1 00252 #define CHOLMOD_VERSION \ 00253 CHOLMOD_VER_CODE(CHOLMOD_MAIN_VERSION,CHOLMOD_SUB_VERSION) 00254 00255 00256 /* ========================================================================== */ 00257 /* === non-CHOLMOD include files ============================================ */ 00258 /* ========================================================================== */ 00259 00260 /* This is the only non-CHOLMOD include file imposed on the user program. 00261 * It required for size_t definition used here. CHOLMOD itself includes other 00262 * ANSI C89 standard #include files, but does not expose them to the user. 00263 * 00264 * CHOLMOD assumes that your C compiler is ANSI C89 compliant. It does not make 00265 * use of ANSI C99 features. 00266 */ 00267 00268 #include <stddef.h> 00269 #include <stdlib.h> 00270 00271 /* ========================================================================== */ 00272 /* === CHOLMOD objects ====================================================== */ 00273 /* ========================================================================== */ 00274 00275 /* Each CHOLMOD object has its own type code. */ 00276 00277 #define CHOLMOD_COMMON 0 00278 #define CHOLMOD_SPARSE 1 00279 #define CHOLMOD_FACTOR 2 00280 #define CHOLMOD_DENSE 3 00281 #define CHOLMOD_TRIPLET 4 00282 00283 /* ========================================================================== */ 00284 /* === CHOLMOD Common ======================================================= */ 00285 /* ========================================================================== */ 00286 00287 /* itype defines the types of integer used: */ 00288 #define CHOLMOD_INT 0 /* all integer arrays are int */ 00289 #define CHOLMOD_INTLONG 1 /* most are int, some are UF_long */ 00290 #define CHOLMOD_LONG 2 /* all integer arrays are UF_long */ 00291 00292 /* The itype of all parameters for all CHOLMOD routines must match. 00293 * FUTURE WORK: CHOLMOD_INTLONG is not yet supported. 00294 */ 00295 00296 /* dtype defines what the numerical type is (double or float): */ 00297 #define CHOLMOD_DOUBLE 0 /* all numerical values are double */ 00298 #define CHOLMOD_SINGLE 1 /* all numerical values are float */ 00299 00300 /* The dtype of all parameters for all CHOLMOD routines must match. 00301 * 00302 * Scalar floating-point values are always passed as double arrays of size 2 00303 * (for the real and imaginary parts). They are typecast to float as needed. 00304 * FUTURE WORK: the float case is not supported yet. 00305 */ 00306 00307 /* xtype defines the kind of numerical values used: */ 00308 #define CHOLMOD_PATTERN 0 /* pattern only, no numerical values */ 00309 #define CHOLMOD_REAL 1 /* a real matrix */ 00310 #define CHOLMOD_COMPLEX 2 /* a complex matrix (ANSI C99 compatible) */ 00311 #define CHOLMOD_ZOMPLEX 3 /* a complex matrix (MATLAB compatible) */ 00312 00313 /* The xtype of all parameters for all CHOLMOD routines must match. 00314 * 00315 * CHOLMOD_PATTERN: x and z are ignored. 00316 * CHOLMOD_DOUBLE: x is non-null of size nzmax, z is ignored. 00317 * CHOLMOD_COMPLEX: x is non-null of size 2*nzmax doubles, z is ignored. 00318 * CHOLMOD_ZOMPLEX: x and z are non-null of size nzmax 00319 * 00320 * In the real case, z is ignored. The kth entry in the matrix is x [k]. 00321 * There are two methods for the complex case. In the ANSI C99-compatible 00322 * CHOLMOD_COMPLEX case, the real and imaginary parts of the kth entry 00323 * are in x [2*k] and x [2*k+1], respectively. z is ignored. In the 00324 * MATLAB-compatible CHOLMOD_ZOMPLEX case, the real and imaginary 00325 * parts of the kth entry are in x [k] and z [k]. 00326 * 00327 * Scalar floating-point values are always passed as double arrays of size 2 00328 * (real and imaginary parts). The imaginary part of a scalar is ignored if 00329 * the routine operates on a real matrix. 00330 * 00331 * These Modules support complex and zomplex matrices, with a few exceptions: 00332 * 00333 * Check all routines 00334 * Cholesky all routines 00335 * Core all except cholmod_aat, add, band, copy 00336 * Demo all routines 00337 * Partition all routines 00338 * Supernodal all routines support any real, complex, or zomplex input. 00339 * There will never be a supernodal zomplex L; a complex 00340 * supernodal L is created if A is zomplex. 00341 * Tcov all routines 00342 * Valgrind all routines 00343 * 00344 * These Modules provide partial support for complex and zomplex matrices: 00345 * 00346 * MATLAB all routines support real and zomplex only, not complex, 00347 * with the exception of ldlupdate, which supports 00348 * real matrices only. This is a minor constraint since 00349 * MATLAB's matrices are all real or zomplex. 00350 * MatrixOps only norm_dense, norm_sparse, and sdmult support complex 00351 * and zomplex 00352 * 00353 * These Modules do not support complex and zomplex matrices at all: 00354 * 00355 * Modify all routines support real matrices only 00356 */ 00357 00358 /* Definitions for cholmod_common: */ 00359 #define CHOLMOD_MAXMETHODS 9 /* maximum number of different methods that */ 00360 /* cholmod_analyze can try. Must be >= 9. */ 00361 00362 /* Common->status values. zero means success, negative means a fatal error, 00363 * positive is a warning. */ 00364 #define CHOLMOD_OK 0 /* success */ 00365 #define CHOLMOD_NOT_INSTALLED (-1) /* failure: method not installed */ 00366 #define CHOLMOD_OUT_OF_MEMORY (-2) /* failure: out of memory */ 00367 #define CHOLMOD_TOO_LARGE (-3) /* failure: integer overflow occured */ 00368 #define CHOLMOD_INVALID (-4) /* failure: invalid input */ 00369 #define CHOLMOD_NOT_POSDEF (1) /* warning: matrix not pos. def. */ 00370 #define CHOLMOD_DSMALL (2) /* warning: D for LDL' or diag(L) or */ 00371 /* LL' has tiny absolute value */ 00372 00373 /* ordering method (also used for L->ordering) */ 00374 #define CHOLMOD_NATURAL 0 /* use natural ordering */ 00375 #define CHOLMOD_GIVEN 1 /* use given permutation */ 00376 #define CHOLMOD_AMD 2 /* use minimum degree (AMD) */ 00377 #define CHOLMOD_METIS 3 /* use METIS' nested dissection */ 00378 #define CHOLMOD_NESDIS 4 /* use CHOLMOD's version of nested dissection:*/ 00379 /* node bisector applied recursively, followed 00380 * by constrained minimum degree (CSYMAMD or 00381 * CCOLAMD) */ 00382 #define CHOLMOD_COLAMD 5 /* use AMD for A, COLAMD for A*A' */ 00383 00384 /* POSTORDERED is not a method, but a result of natural ordering followed by a 00385 * weighted postorder. It is used for L->ordering, not method [ ].ordering. */ 00386 #define CHOLMOD_POSTORDERED 6 /* natural ordering, postordered. */ 00387 00388 /* supernodal strategy (for Common->supernodal) */ 00389 #define CHOLMOD_SIMPLICIAL 0 /* always do simplicial */ 00390 #define CHOLMOD_AUTO 1 /* select simpl/super depending on matrix */ 00391 #define CHOLMOD_SUPERNODAL 2 /* always do supernodal */ 00392 00393 typedef struct cholmod_common_struct 00394 { 00395 /* ---------------------------------------------------------------------- */ 00396 /* parameters for symbolic/numeric factorization and update/downdate */ 00397 /* ---------------------------------------------------------------------- */ 00398 00399 double dbound ; /* Smallest absolute value of diagonal entries of D 00400 * for LDL' factorization and update/downdate/rowadd/ 00401 * rowdel, or the diagonal of L for an LL' factorization. 00402 * Entries in the range 0 to dbound are replaced with dbound. 00403 * Entries in the range -dbound to 0 are replaced with -dbound. No 00404 * changes are made to the diagonal if dbound <= 0. Default: zero */ 00405 00406 double grow0 ; /* For a simplicial factorization, L->i and L->x can 00407 * grow if necessary. grow0 is the factor by which 00408 * it grows. For the initial space, L is of size MAX (1,grow0) times 00409 * the required space. If L runs out of space, the new size of L is 00410 * MAX(1.2,grow0) times the new required space. If you do not plan on 00411 * modifying the LDL' factorization in the Modify module, set grow0 to 00412 * zero (or set grow2 to 0, see below). Default: 1.2 */ 00413 00414 double grow1 ; 00415 00416 size_t grow2 ; /* For a simplicial factorization, each column j of L 00417 * is initialized with space equal to 00418 * grow1*L->ColCount[j] + grow2. If grow0 < 1, grow1 < 1, or grow2 == 0, 00419 * then the space allocated is exactly equal to L->ColCount[j]. If the 00420 * column j runs out of space, it increases to grow1*need + grow2 in 00421 * size, where need is the total # of nonzeros in that column. If you do 00422 * not plan on modifying the factorization in the Modify module, set 00423 * grow2 to zero. Default: grow1 = 1.2, grow2 = 5. */ 00424 00425 size_t maxrank ; /* rank of maximum update/downdate. Valid values: 00426 * 2, 4, or 8. A value < 2 is set to 2, and a 00427 * value > 8 is set to 8. It is then rounded up to the next highest 00428 * power of 2, if not already a power of 2. Workspace (Xwork, below) of 00429 * size nrow-by-maxrank double's is allocated for the update/downdate. 00430 * If an update/downdate of rank-k is requested, with k > maxrank, 00431 * it is done in steps of maxrank. Default: 8, which is fastest. 00432 * Memory usage can be reduced by setting maxrank to 2 or 4. 00433 */ 00434 00435 double supernodal_switch ; /* supernodal vs simplicial factorization */ 00436 int supernodal ; /* If Common->supernodal <= CHOLMOD_SIMPLICIAL 00437 * (0) then cholmod_analyze performs a 00438 * simplicial analysis. If >= CHOLMOD_SUPERNODAL (2), then a supernodal 00439 * analysis is performed. If == CHOLMOD_AUTO (1) and 00440 * flop/nnz(L) < Common->supernodal_switch, then a simplicial analysis 00441 * is done. A supernodal analysis done otherwise. 00442 * Default: CHOLMOD_AUTO. Default supernodal_switch = 40 */ 00443 00444 int final_asis ; /* If TRUE, then ignore the other final_* parameters 00445 * (except for final_pack). 00446 * The factor is left as-is when done. Default: TRUE.*/ 00447 00448 int final_super ; /* If TRUE, leave a factor in supernodal form when 00449 * supernodal factorization is finished. If FALSE, 00450 * then convert to a simplicial factor when done. 00451 * Default: TRUE */ 00452 00453 int final_ll ; /* If TRUE, leave factor in LL' form when done. 00454 * Otherwise, leave in LDL' form. Default: FALSE */ 00455 00456 int final_pack ; /* If TRUE, pack the columns when done. If TRUE, and 00457 * cholmod_factorize is called with a symbolic L, L is 00458 * allocated with exactly the space required, using L->ColCount. If you 00459 * plan on modifying the factorization, set Common->final_pack to FALSE, 00460 * and each column will be given a little extra slack space for future 00461 * growth in fill-in due to updates. Default: TRUE */ 00462 00463 int final_monotonic ; /* If TRUE, ensure columns are monotonic when done. 00464 * Default: TRUE */ 00465 00466 int final_resymbol ;/* if cholmod_factorize performed a supernodal 00467 * factorization, final_resymbol is true, and 00468 * final_super is FALSE (convert a simplicial numeric factorization), 00469 * then numerically zero entries that resulted from relaxed supernodal 00470 * amalgamation are removed. This does not remove entries that are zero 00471 * due to exact numeric cancellation, since doing so would break the 00472 * update/downdate rowadd/rowdel routines. Default: FALSE. */ 00473 00474 /* supernodal relaxed amalgamation parameters: */ 00475 double zrelax [3] ; 00476 size_t nrelax [3] ; 00477 00478 /* Let ns be the total number of columns in two adjacent supernodes. 00479 * Let z be the fraction of zero entries in the two supernodes if they 00480 * are merged (z includes zero entries from prior amalgamations). The 00481 * two supernodes are merged if: 00482 * (ns <= nrelax [0]) || (no new zero entries added) || 00483 * (ns <= nrelax [1] && z < zrelax [0]) || 00484 * (ns <= nrelax [2] && z < zrelax [1]) || (z < zrelax [2]) 00485 * 00486 * Default parameters result in the following rule: 00487 * (ns <= 4) || (no new zero entries added) || 00488 * (ns <= 16 && z < 0.8) || (ns <= 48 && z < 0.1) || (z < 0.05) 00489 */ 00490 00491 int prefer_zomplex ; /* X = cholmod_solve (sys, L, B, Common) computes 00492 * x=A\b or solves a related system. If L and B are 00493 * both real, then X is real. Otherwise, X is returned as 00494 * CHOLMOD_COMPLEX if Common->prefer_zomplex is FALSE, or 00495 * CHOLMOD_ZOMPLEX if Common->prefer_zomplex is TRUE. This parameter 00496 * is needed because there is no supernodal zomplex L. Suppose the 00497 * caller wants all complex matrices to be stored in zomplex form 00498 * (MATLAB, for example). A supernodal L is returned in complex form 00499 * if A is zomplex. B can be real, and thus X = cholmod_solve (L,B) 00500 * should return X as zomplex. This cannot be inferred from the input 00501 * arguments L and B. Default: FALSE, since all data types are 00502 * supported in CHOLMOD_COMPLEX form and since this is the native type 00503 * of LAPACK and the BLAS. Note that the MATLAB/cholmod.c mexFunction 00504 * sets this parameter to TRUE, since MATLAB matrices are in 00505 * CHOLMOD_ZOMPLEX form. 00506 */ 00507 00508 int prefer_upper ; /* cholmod_analyze and cholmod_factorize work 00509 * fastest when a symmetric matrix is stored in 00510 * upper triangular form when a fill-reducing ordering is used. In 00511 * MATLAB, this corresponds to how x=A\b works. When the matrix is 00512 * ordered as-is, they work fastest when a symmetric matrix is in lower 00513 * triangular form. In MATLAB, R=chol(A) does the opposite. This 00514 * parameter affects only how cholmod_read returns a symmetric matrix. 00515 * If TRUE (the default case), a symmetric matrix is always returned in 00516 * upper-triangular form (A->stype = 1). */ 00517 00518 int quick_return_if_not_posdef ; /* if TRUE, the supernodal numeric 00519 * factorization will return quickly if 00520 * the matrix is not positive definite. Default: FALSE. */ 00521 00522 /* ---------------------------------------------------------------------- */ 00523 /* printing and error handling options */ 00524 /* ---------------------------------------------------------------------- */ 00525 00526 int print ; /* print level. Default: 3 */ 00527 int precise ; /* if TRUE, print 16 digits. Otherwise print 5 */ 00528 int (*print_function) (const char *, ...) ; /* pointer to printf */ 00529 00530 int try_catch ; /* if TRUE, then ignore errors; CHOLMOD is in the middle 00531 * of a try/catch block. No error message is printed 00532 * and the Common->error_handler function is not called. */ 00533 00534 void (*error_handler) (int status, const char *file, 00535 int line, const char *message) ; 00536 00537 /* Common->error_handler is the user's error handling routine. If not 00538 * NULL, this routine is called if an error occurs in CHOLMOD. status 00539 * can be CHOLMOD_OK (0), negative for a fatal error, and positive for 00540 * a warning. file is a string containing the name of the source code 00541 * file where the error occured, and line is the line number in that 00542 * file. message is a string describing the error in more detail. */ 00543 00544 /* ---------------------------------------------------------------------- */ 00545 /* ordering options */ 00546 /* ---------------------------------------------------------------------- */ 00547 00548 /* The cholmod_analyze routine can try many different orderings and select 00549 * the best one. It can also try one ordering method multiple times, with 00550 * different parameter settings. The default is to use three orderings, 00551 * the user's permutation (if provided), AMD which is the fastest ordering 00552 * and generally gives good fill-in, and METIS. CHOLMOD's nested dissection 00553 * (METIS with a constrained AMD) usually gives a better ordering than METIS 00554 * alone (by about 5% to 10%) but it takes more time. 00555 * 00556 * If you know the method that is best for your matrix, set Common->nmethods 00557 * to 1 and set Common->method [0] to the set of parameters for that method. 00558 * If you set it to 1 and do not provide a permutation, then only AMD will 00559 * be called. 00560 * 00561 * If METIS is not available, the default # of methods tried is 2 (the user 00562 * permutation, if any, and AMD). 00563 * 00564 * To try other methods, set Common->nmethods to the number of methods you 00565 * want to try. The suite of default methods and their parameters is 00566 * described in the cholmod_defaults routine, and summarized here: 00567 * 00568 * Common->method [i]: 00569 * i = 0: user-provided ordering (cholmod_analyze_p only) 00570 * i = 1: AMD (for both A and A*A') 00571 * i = 2: METIS 00572 * i = 3: CHOLMOD's nested dissection (NESDIS), default parameters 00573 * i = 4: natural 00574 * i = 5: NESDIS with nd_small = 20000 00575 * i = 6: NESDIS with nd_small = 4, no constrained minimum degree 00576 * i = 7: NESDIS with no dense node removal 00577 * i = 8: AMD for A, COLAMD for A*A' 00578 * 00579 * You can modify the suite of methods you wish to try by modifying 00580 * Common.method [...] after calling cholmod_start or cholmod_defaults. 00581 * 00582 * For example, to use AMD, followed by a weighted postordering: 00583 * 00584 * Common->nmethods = 1 ; 00585 * Common->method [0].ordering = CHOLMOD_AMD ; 00586 * Common->postorder = TRUE ; 00587 * 00588 * To use the natural ordering (with no postordering): 00589 * 00590 * Common->nmethods = 1 ; 00591 * Common->method [0].ordering = CHOLMOD_NATURAL ; 00592 * Common->postorder = FALSE ; 00593 * 00594 * If you are going to factorize hundreds or more matrices with the same 00595 * nonzero pattern, you may wish to spend a great deal of time finding a 00596 * good permutation. In this case, try setting Common->nmethods to 9. 00597 * The time spent in cholmod_analysis will be very high, but you need to 00598 * call it only once. 00599 * 00600 * cholmod_analyze sets Common->current to a value between 0 and nmethods-1. 00601 * Each ordering method uses the set of options defined by this parameter. 00602 */ 00603 00604 int nmethods ; /* The number of ordering methods to try. Default: 0. 00605 * nmethods = 0 is a special case. cholmod_analyze 00606 * will try the user-provided ordering (if given) and AMD. Let fl and 00607 * lnz be the flop count and nonzeros in L from AMD's ordering. Let 00608 * anz be the number of nonzeros in the upper or lower triangular part 00609 * of the symmetric matrix A. If fl/lnz < 500 or lnz/anz < 5, then this 00610 * is a good ordering, and METIS is not attempted. Otherwise, METIS is 00611 * tried. The best ordering found is used. If nmethods > 0, the 00612 * methods used are given in the method[ ] array, below. The first 00613 * three methods in the default suite of orderings is (1) use the given 00614 * permutation (if provided), (2) use AMD, and (3) use METIS. Maximum 00615 * allowed value is CHOLMOD_MAXMETHODS. */ 00616 00617 int current ; /* The current method being tried. Default: 0. Valid 00618 * range is 0 to nmethods-1. */ 00619 00620 int selected ; /* The best method found. */ 00621 00622 /* The suite of ordering methods and parameters: */ 00623 00624 struct cholmod_method_struct 00625 { 00626 /* statistics for this method */ 00627 double lnz ; /* nnz(L) excl. zeros from supernodal amalgamation, 00628 * for a "pure" L */ 00629 00630 double fl ; /* flop count for a "pure", real simplicial LL' 00631 * factorization, with no extra work due to 00632 * amalgamation. Subtract n to get the LDL' flop count. Multiply 00633 * by about 4 if the matrix is complex or zomplex. */ 00634 00635 /* ordering method parameters */ 00636 double prune_dense ;/* dense row/col control for AMD, SYMAMD, CSYMAMD, 00637 * and NESDIS (cholmod_nested_dissection). For a 00638 * symmetric n-by-n matrix, rows/columns with more than 00639 * MAX (16, prune_dense * sqrt (n)) entries are removed prior to 00640 * ordering. They appear at the end of the re-ordered matrix. 00641 * 00642 * If prune_dense < 0, only completely dense rows/cols are removed. 00643 * 00644 * This paramater is also the dense column control for COLAMD and 00645 * CCOLAMD. For an m-by-n matrix, columns with more than 00646 * MAX (16, prune_dense * sqrt (MIN (m,n))) entries are removed prior 00647 * to ordering. They appear at the end of the re-ordered matrix. 00648 * CHOLMOD factorizes A*A', so it calls COLAMD and CCOLAMD with A', 00649 * not A. Thus, this parameter affects the dense *row* control for 00650 * CHOLMOD's matrix, and the dense *column* control for COLAMD and 00651 * CCOLAMD. 00652 * 00653 * Removing dense rows and columns improves the run-time of the 00654 * ordering methods. It has some impact on ordering quality 00655 * (usually minimal, sometimes good, sometimes bad). 00656 * 00657 * Default: 10. */ 00658 00659 double prune_dense2 ;/* dense row control for COLAMD and CCOLAMD. 00660 * Rows with more than MAX (16, dense2 * sqrt (n)) 00661 * for an m-by-n matrix are removed prior to ordering. CHOLMOD's 00662 * matrix is transposed before ordering it with COLAMD or CCOLAMD, 00663 * so this controls the dense *columns* of CHOLMOD's matrix, and 00664 * the dense *rows* of COLAMD's or CCOLAMD's matrix. 00665 * 00666 * If prune_dense2 < 0, only completely dense rows/cols are removed. 00667 * 00668 * Default: -1. Note that this is not the default for COLAMD and 00669 * CCOLAMD. -1 is best for Cholesky. 10 is best for LU. */ 00670 00671 double nd_oksep ; /* in NESDIS, when a node separator is computed, it 00672 * discarded if nsep >= nd_oksep*n, where nsep is 00673 * the number of nodes in the separator, and n is the size of the 00674 * graph being cut. Valid range is 0 to 1. If 1 or greater, the 00675 * separator is discarded if it consists of the entire graph. 00676 * Default: 1 */ 00677 00678 double other1 [4] ; /* future expansion */ 00679 00680 size_t nd_small ; /* do not partition graphs with fewer nodes than 00681 * nd_small, in NESDIS. Default: 200 (same as 00682 * METIS) */ 00683 00684 size_t other2 [4] ; /* future expansion */ 00685 00686 int aggressive ; /* Aggresive absorption in AMD, COLAMD, SYMAMD, 00687 * CCOLAMD, and CSYMAMD. Default: TRUE */ 00688 00689 int order_for_lu ; /* CCOLAMD can be optimized to produce an ordering 00690 * for LU or Cholesky factorization. CHOLMOD only 00691 * performs a Cholesky factorization. However, you may wish to use 00692 * CHOLMOD as an interface for CCOLAMD but use it for your own LU 00693 * factorization. In this case, order_for_lu should be set to FALSE. 00694 * When factorizing in CHOLMOD itself, you should *** NEVER *** set 00695 * this parameter FALSE. Default: TRUE. */ 00696 00697 int nd_compress ; /* If TRUE, compress the graph and subgraphs before 00698 * partitioning them in NESDIS. Default: TRUE */ 00699 00700 int nd_camd ; /* If 1, follow the nested dissection ordering 00701 * with a constrained minimum degree ordering that 00702 * respects the partitioning just found (using CAMD). If 2, use 00703 * CSYMAMD instead. If you set nd_small very small, you may not need 00704 * this ordering, and can save time by setting it to zero (no 00705 * constrained minimum degree ordering). Default: 1. */ 00706 00707 int nd_components ; /* The nested dissection ordering finds a node 00708 * separator that splits the graph into two parts, 00709 * which may be unconnected. If nd_components is TRUE, each of 00710 * these connected components is split independently. If FALSE, 00711 * each part is split as a whole, even if it consists of more than 00712 * one connected component. Default: FALSE */ 00713 00714 /* fill-reducing ordering to use */ 00715 int ordering ; 00716 00717 size_t other3 [4] ; /* future expansion */ 00718 00719 } method [CHOLMOD_MAXMETHODS + 1] ; 00720 00721 int postorder ; /* If TRUE, cholmod_analyze follows the ordering with a 00722 * weighted postorder of the elimination tree. Improves 00723 * supernode amalgamation. Does not affect fundamental nnz(L) and 00724 * flop count. Default: TRUE. */ 00725 00726 /* ---------------------------------------------------------------------- */ 00727 /* memory management routines */ 00728 /* ---------------------------------------------------------------------- */ 00729 00730 void *(*malloc_memory) (size_t) ; /* pointer to malloc */ 00731 void *(*realloc_memory) (void *, size_t) ; /* pointer to realloc */ 00732 void (*free_memory) (void *) ; /* pointer to free */ 00733 void *(*calloc_memory) (size_t, size_t) ; /* pointer to calloc */ 00734 00735 /* ---------------------------------------------------------------------- */ 00736 /* routines for complex arithmetic */ 00737 /* ---------------------------------------------------------------------- */ 00738 00739 int (*complex_divide) (double ax, double az, double bx, double bz, 00740 double *cx, double *cz) ; 00741 00742 /* flag = complex_divide (ax, az, bx, bz, &cx, &cz) computes the complex 00743 * division c = a/b, where ax and az hold the real and imaginary part 00744 * of a, and b and c are stored similarly. flag is returned as 1 if 00745 * a divide-by-zero occurs, or 0 otherwise. By default, the function 00746 * pointer Common->complex_divide is set equal to cholmod_divcomplex. 00747 */ 00748 00749 double (*hypotenuse) (double x, double y) ; 00750 00751 /* s = hypotenuse (x,y) computes s = sqrt (x*x + y*y), but does so more 00752 * accurately. By default, the function pointer Common->hypotenuse is 00753 * set equal to cholmod_hypot. See also the hypot function in the C99 00754 * standard, which has an identical syntax and function. If you have 00755 * a C99-compliant compiler, you can set Common->hypotenuse = hypot. */ 00756 00757 /* ---------------------------------------------------------------------- */ 00758 /* METIS workarounds */ 00759 /* ---------------------------------------------------------------------- */ 00760 00761 double metis_memory ; /* This is a parameter for CHOLMOD's interface to 00762 * METIS, not a parameter to METIS itself. METIS 00763 * uses an amount of memory that is difficult to estimate precisely 00764 * beforehand. If it runs out of memory, it terminates your program. 00765 * All routines in CHOLMOD except for CHOLMOD's interface to METIS 00766 * return an error status and safely return to your program if they run 00767 * out of memory. To mitigate this problem, the CHOLMOD interface 00768 * can allocate a single block of memory equal in size to an empirical 00769 * upper bound of METIS's memory usage times the Common->metis_memory 00770 * parameter, and then immediately free it. It then calls METIS. If 00771 * this pre-allocation fails, it is possible that METIS will fail as 00772 * well, and so CHOLMOD returns with an out-of-memory condition without 00773 * calling METIS. 00774 * 00775 * METIS_NodeND (used in the CHOLMOD_METIS ordering option) with its 00776 * default parameter settings typically uses about (4*nz+40n+4096) 00777 * times sizeof(int) memory, where nz is equal to the number of entries 00778 * in A for the symmetric case or AA' if an unsymmetric matrix is 00779 * being ordered (where nz includes both the upper and lower parts 00780 * of A or AA'). The observed "upper bound" (with 2 exceptions), 00781 * measured in an instrumented copy of METIS 4.0.1 on thousands of 00782 * matrices, is (10*nz+50*n+4096) * sizeof(int). Two large matrices 00783 * exceeded this bound, one by almost a factor of 2 (Gupta/gupta2). 00784 * 00785 * If your program is terminated by METIS, try setting metis_memory to 00786 * 2.0, or even higher if needed. By default, CHOLMOD assumes that METIS 00787 * does not have this problem (so that CHOLMOD will work correctly when 00788 * this issue is fixed in METIS). Thus, the default value is zero. 00789 * This work-around is not guaranteed anyway. 00790 * 00791 * If a matrix exceeds this predicted memory usage, AMD is attempted 00792 * instead. It, too, may run out of memory, but if it does so it will 00793 * not terminate your program. 00794 */ 00795 00796 double metis_dswitch ; /* METIS_NodeND in METIS 4.0.1 gives a seg */ 00797 size_t metis_nswitch ; /* fault with one matrix of order n = 3005 and 00798 * nz = 6,036,025. This is a very dense graph. 00799 * The workaround is to use AMD instead of METIS for matrices of dimension 00800 * greater than Common->metis_nswitch (default 3000) or more and with 00801 * density of Common->metis_dswitch (default 0.66) or more. 00802 * cholmod_nested_dissection has no problems with the same matrix, even 00803 * though it uses METIS_NodeComputeSeparator on this matrix. If this 00804 * seg fault does not affect you, set metis_nswitch to zero or less, 00805 * and CHOLMOD will not switch to AMD based just on the density of the 00806 * matrix (it will still switch to AMD if the metis_memory parameter 00807 * causes the switch). 00808 */ 00809 00810 /* ---------------------------------------------------------------------- */ 00811 /* workspace */ 00812 /* ---------------------------------------------------------------------- */ 00813 00814 /* CHOLMOD has several routines that take less time than the size of 00815 * workspace they require. Allocating and initializing the workspace would 00816 * dominate the run time, unless workspace is allocated and initialized 00817 * just once. CHOLMOD allocates this space when needed, and holds it here 00818 * between calls to CHOLMOD. cholmod_start sets these pointers to NULL 00819 * (which is why it must be the first routine called in CHOLMOD). 00820 * cholmod_finish frees the workspace (which is why it must be the last 00821 * call to CHOLMOD). 00822 */ 00823 00824 size_t nrow ; /* size of Flag and Head */ 00825 UF_long mark ; /* mark value for Flag array */ 00826 size_t iworksize ; /* size of Iwork. Upper bound: 6*nrow+ncol */ 00827 size_t xworksize ; /* size of Xwork, in bytes. 00828 * maxrank*nrow*sizeof(double) for update/downdate. 00829 * 2*nrow*sizeof(double) otherwise */ 00830 00831 /* initialized workspace: contents needed between calls to CHOLMOD */ 00832 void *Flag ; /* size nrow, an integer array. Kept cleared between 00833 * calls to cholmod rouines (Flag [i] < mark) */ 00834 00835 void *Head ; /* size nrow+1, an integer array. Kept cleared between 00836 * calls to cholmod routines (Head [i] = EMPTY) */ 00837 00838 void *Xwork ; /* a double array. Its size varies. It is nrow for 00839 * most routines (cholmod_rowfac, cholmod_add, 00840 * cholmod_aat, cholmod_norm, cholmod_ssmult) for the real case, twice 00841 * that when the input matrices are complex or zomplex. It is of size 00842 * 2*nrow for cholmod_rowadd and cholmod_rowdel. For cholmod_updown, 00843 * its size is maxrank*nrow where maxrank is 2, 4, or 8. Kept cleared 00844 * between calls to cholmod (set to zero). */ 00845 00846 /* uninitialized workspace, contents not needed between calls to CHOLMOD */ 00847 void *Iwork ; /* size iworksize, 2*nrow+ncol for most routines, 00848 * up to 6*nrow+ncol for cholmod_analyze. */ 00849 00850 int itype ; /* If CHOLMOD_LONG, Flag, Head, and Iwork are UF_long. 00851 * Otherwise all three arrays are int. */ 00852 00853 int dtype ; /* double or float */ 00854 00855 /* Common->itype and Common->dtype are used to define the types of all 00856 * sparse matrices, triplet matrices, dense matrices, and factors 00857 * created using this Common struct. The itypes and dtypes of all 00858 * parameters to all CHOLMOD routines must match. */ 00859 00860 int no_workspace_reallocate ; /* this is an internal flag, used as a 00861 * precaution by cholmod_analyze. It is normally false. If true, 00862 * cholmod_allocate_work is not allowed to reallocate any workspace; 00863 * they must use the existing workspace in Common (Iwork, Flag, Head, 00864 * and Xwork). Added for CHOLMOD v1.1 */ 00865 00866 /* ---------------------------------------------------------------------- */ 00867 /* statistics */ 00868 /* ---------------------------------------------------------------------- */ 00869 00870 /* fl and lnz are set only in cholmod_analyze and cholmod_rowcolcounts, 00871 * in the Cholesky modudle. modfl is set only in the Modify module. */ 00872 00873 int status ; /* error code */ 00874 double fl ; /* LL' flop count from most recent analysis */ 00875 double lnz ; /* fundamental nz in L */ 00876 double anz ; /* nonzeros in tril(A) if A is symmetric/lower, 00877 * triu(A) if symmetric/upper, or tril(A*A') if 00878 * unsymmetric, in last call to cholmod_analyze. */ 00879 double modfl ; /* flop count from most recent update/downdate/ 00880 * rowadd/rowdel (excluding flops to modify the 00881 * solution to Lx=b, if computed) */ 00882 size_t malloc_count ; /* # of objects malloc'ed minus the # free'd*/ 00883 size_t memory_usage ; /* peak memory usage in bytes */ 00884 size_t memory_inuse ; /* current memory usage in bytes */ 00885 00886 double nrealloc_col ; /* # of column reallocations */ 00887 double nrealloc_factor ;/* # of factor reallocations due to col. reallocs */ 00888 double ndbounds_hit ; /* # of times diagonal modified by dbound */ 00889 00890 double rowfacfl ; /* # of flops in last call to cholmod_rowfac */ 00891 double aatfl ; /* # of flops to compute A(:,f)*A(:,f)' */ 00892 00893 /* ---------------------------------------------------------------------- */ 00894 /* future expansion */ 00895 /* ---------------------------------------------------------------------- */ 00896 00897 /* To allow CHOLMOD to be updated without recompiling the user application, 00898 * additional space is set aside here for future statistics, parameters, 00899 * and workspace. Note: additional entries were added in v1.1 to the 00900 * method array, above, and thus v1.0 and v1.1 are not binary compatible. 00901 * 00902 * v1.1 to the current version are binary compatible. 00903 */ 00904 00905 /* ---------------------------------------------------------------------- */ 00906 double other1 [12] ; /* reduced from size 16 in v1.6 */ 00907 00908 double SPQR_xstat [2] ; /* for SuiteSparseQR statistics */ 00909 00910 /* SuiteSparseQR control parameters: */ 00911 double SPQR_grain ; /* task size is >= max (total flops / grain) */ 00912 double SPQR_small ; /* task size is >= small */ 00913 00914 /* ---------------------------------------------------------------------- */ 00915 UF_long SPQR_istat [10] ; /* for SuiteSparseQR statistics */ 00916 UF_long other2 [6] ; /* reduced from size 16 in v1.6 */ 00917 00918 /* ---------------------------------------------------------------------- */ 00919 int other3 [10] ; /* reduced from size 16 in v1.1. */ 00920 00921 int prefer_binary ; /* cholmod_read_triplet converts a symmetric 00922 * pattern-only matrix into a real matrix. If 00923 * prefer_binary is FALSE, the diagonal entries are set to 1 + the degree 00924 * of the row/column, and off-diagonal entries are set to -1 (resulting 00925 * in a positive definite matrix if the diagonal is zero-free). Most 00926 * symmetric patterns are the pattern a positive definite matrix. If 00927 * this parameter is TRUE, then the matrix is returned with a 1 in each 00928 * entry, instead. Default: FALSE. Added in v1.3. */ 00929 00930 /* control parameter (added for v1.2): */ 00931 int default_nesdis ; /* Default: FALSE. If FALSE, then the default 00932 * ordering strategy (when Common->nmethods == 0) 00933 * is to try the given ordering (if present), AMD, and then METIS if AMD 00934 * reports high fill-in. If Common->default_nesdis is TRUE then NESDIS 00935 * is used instead in the default strategy. */ 00936 00937 /* statistic (added for v1.2): */ 00938 int called_nd ; /* TRUE if the last call to 00939 * cholmod_analyze called NESDIS or METIS. */ 00940 00941 int blas_ok ; /* FALSE if BLAS int overflow; TRUE otherwise */ 00942 00943 /* SuiteSparseQR control parameters: */ 00944 int SPQR_shrink ; /* controls stack realloc method */ 00945 int SPQR_nthreads ; /* number of TBB threads, 0 = auto */ 00946 00947 /* ---------------------------------------------------------------------- */ 00948 size_t other4 [16] ; 00949 00950 /* ---------------------------------------------------------------------- */ 00951 void *other5 [16] ; 00952 00953 } cholmod_common ; 00954 00955 /* -------------------------------------------------------------------------- */ 00956 /* cholmod_start: first call to CHOLMOD */ 00957 /* -------------------------------------------------------------------------- */ 00958 00959 int cholmod_start 00960 ( 00961 cholmod_common *Common 00962 ) ; 00963 00964 int cholmod_l_start (cholmod_common *) ; 00965 00966 /* -------------------------------------------------------------------------- */ 00967 /* cholmod_finish: last call to CHOLMOD */ 00968 /* -------------------------------------------------------------------------- */ 00969 00970 int cholmod_finish 00971 ( 00972 cholmod_common *Common 00973 ) ; 00974 00975 int cholmod_l_finish (cholmod_common *) ; 00976 00977 /* -------------------------------------------------------------------------- */ 00978 /* cholmod_defaults: restore default parameters */ 00979 /* -------------------------------------------------------------------------- */ 00980 00981 int cholmod_defaults 00982 ( 00983 cholmod_common *Common 00984 ) ; 00985 00986 int cholmod_l_defaults (cholmod_common *) ; 00987 00988 /* -------------------------------------------------------------------------- */ 00989 /* cholmod_maxrank: return valid maximum rank for update/downdate */ 00990 /* -------------------------------------------------------------------------- */ 00991 00992 size_t cholmod_maxrank /* returns validated value of Common->maxrank */ 00993 ( 00994 /* ---- input ---- */ 00995 size_t n, /* A and L will have n rows */ 00996 /* --------------- */ 00997 cholmod_common *Common 00998 ) ; 00999 01000 size_t cholmod_l_maxrank (size_t, cholmod_common *) ; 01001 01002 /* -------------------------------------------------------------------------- */ 01003 /* cholmod_allocate_work: allocate workspace in Common */ 01004 /* -------------------------------------------------------------------------- */ 01005 01006 int cholmod_allocate_work 01007 ( 01008 /* ---- input ---- */ 01009 size_t nrow, /* size: Common->Flag (nrow), Common->Head (nrow+1) */ 01010 size_t iworksize, /* size of Common->Iwork */ 01011 size_t xworksize, /* size of Common->Xwork */ 01012 /* --------------- */ 01013 cholmod_common *Common 01014 ) ; 01015 01016 int cholmod_l_allocate_work (size_t, size_t, size_t, cholmod_common *) ; 01017 01018 /* -------------------------------------------------------------------------- */ 01019 /* cholmod_free_work: free workspace in Common */ 01020 /* -------------------------------------------------------------------------- */ 01021 01022 int cholmod_free_work 01023 ( 01024 cholmod_common *Common 01025 ) ; 01026 01027 int cholmod_l_free_work (cholmod_common *) ; 01028 01029 /* -------------------------------------------------------------------------- */ 01030 /* cholmod_clear_flag: clear Flag workspace in Common */ 01031 /* -------------------------------------------------------------------------- */ 01032 01033 /* use a macro for speed */ 01034 #define CHOLMOD_CLEAR_FLAG(Common) \ 01035 { \ 01036 Common->mark++ ; \ 01037 if (Common->mark <= 0) \ 01038 { \ 01039 Common->mark = EMPTY ; \ 01040 CHOLMOD (clear_flag) (Common) ; \ 01041 } \ 01042 } 01043 01044 UF_long cholmod_clear_flag 01045 ( 01046 cholmod_common *Common 01047 ) ; 01048 01049 UF_long cholmod_l_clear_flag (cholmod_common *) ; 01050 01051 /* -------------------------------------------------------------------------- */ 01052 /* cholmod_error: called when CHOLMOD encounters an error */ 01053 /* -------------------------------------------------------------------------- */ 01054 01055 int cholmod_error 01056 ( 01057 /* ---- input ---- */ 01058 int status, /* error status */ 01059 const char *file, /* name of source code file where error occured */ 01060 int line, /* line number in source code file where error occured*/ 01061 const char *message,/* error message */ 01062 /* --------------- */ 01063 cholmod_common *Common 01064 ) ; 01065 01066 int cholmod_l_error (int, const char *, int, const char *, cholmod_common *) ; 01067 01068 /* -------------------------------------------------------------------------- */ 01069 /* cholmod_dbound: for internal use in CHOLMOD only */ 01070 /* -------------------------------------------------------------------------- */ 01071 01072 double cholmod_dbound /* returns modified diagonal entry of D or L */ 01073 ( 01074 /* ---- input ---- */ 01075 double dj, /* diagonal entry of D for LDL' or L for LL' */ 01076 /* --------------- */ 01077 cholmod_common *Common 01078 ) ; 01079 01080 double cholmod_l_dbound (double, cholmod_common *) ; 01081 01082 /* -------------------------------------------------------------------------- */ 01083 /* cholmod_hypot: compute sqrt (x*x + y*y) accurately */ 01084 /* -------------------------------------------------------------------------- */ 01085 01086 double cholmod_hypot 01087 ( 01088 /* ---- input ---- */ 01089 double x, double y 01090 ) ; 01091 01092 double cholmod_l_hypot (double, double) ; 01093 01094 /* -------------------------------------------------------------------------- */ 01095 /* cholmod_divcomplex: complex division, c = a/b */ 01096 /* -------------------------------------------------------------------------- */ 01097 01098 int cholmod_divcomplex /* return 1 if divide-by-zero, 0 otherise */ 01099 ( 01100 /* ---- input ---- */ 01101 double ar, double ai, /* real and imaginary parts of a */ 01102 double br, double bi, /* real and imaginary parts of b */ 01103 /* ---- output --- */ 01104 double *cr, double *ci /* real and imaginary parts of c */ 01105 ) ; 01106 01107 int cholmod_l_divcomplex (double, double, double, double, double *, double *) ; 01108 01109 01110 /* ========================================================================== */ 01111 /* === Core/cholmod_sparse ================================================== */ 01112 /* ========================================================================== */ 01113 01114 /* A sparse matrix stored in compressed-column form. */ 01115 01116 typedef struct cholmod_sparse_struct 01117 { 01118 size_t nrow ; /* the matrix is nrow-by-ncol */ 01119 size_t ncol ; 01120 size_t nzmax ; /* maximum number of entries in the matrix */ 01121 01122 /* pointers to int or UF_long: */ 01123 void *p ; /* p [0..ncol], the column pointers */ 01124 void *i ; /* i [0..nzmax-1], the row indices */ 01125 01126 /* for unpacked matrices only: */ 01127 void *nz ; /* nz [0..ncol-1], the # of nonzeros in each col. In 01128 * packed form, the nonzero pattern of column j is in 01129 * A->i [A->p [j] ... A->p [j+1]-1]. In unpacked form, column j is in 01130 * A->i [A->p [j] ... A->p [j]+A->nz[j]-1] instead. In both cases, the 01131 * numerical values (if present) are in the corresponding locations in 01132 * the array x (or z if A->xtype is CHOLMOD_ZOMPLEX). */ 01133 01134 /* pointers to double or float: */ 01135 void *x ; /* size nzmax or 2*nzmax, if present */ 01136 void *z ; /* size nzmax, if present */ 01137 01138 int stype ; /* Describes what parts of the matrix are considered: 01139 * 01140 * 0: matrix is "unsymmetric": use both upper and lower triangular parts 01141 * (the matrix may actually be symmetric in pattern and value, but 01142 * both parts are explicitly stored and used). May be square or 01143 * rectangular. 01144 * >0: matrix is square and symmetric, use upper triangular part. 01145 * Entries in the lower triangular part are ignored. 01146 * <0: matrix is square and symmetric, use lower triangular part. 01147 * Entries in the upper triangular part are ignored. 01148 * 01149 * Note that stype>0 and stype<0 are different for cholmod_sparse and 01150 * cholmod_triplet. See the cholmod_triplet data structure for more 01151 * details. 01152 */ 01153 01154 int itype ; /* CHOLMOD_INT: p, i, and nz are int. 01155 * CHOLMOD_INTLONG: p is UF_long, i and nz are int. 01156 * CHOLMOD_LONG: p, i, and nz are UF_long. */ 01157 01158 int xtype ; /* pattern, real, complex, or zomplex */ 01159 int dtype ; /* x and z are double or float */ 01160 int sorted ; /* TRUE if columns are sorted, FALSE otherwise */ 01161 int packed ; /* TRUE if packed (nz ignored), FALSE if unpacked 01162 * (nz is required) */ 01163 01164 } cholmod_sparse ; 01165 01166 /* -------------------------------------------------------------------------- */ 01167 /* cholmod_allocate_sparse: allocate a sparse matrix */ 01168 /* -------------------------------------------------------------------------- */ 01169 01170 cholmod_sparse *cholmod_allocate_sparse 01171 ( 01172 /* ---- input ---- */ 01173 size_t nrow, /* # of rows of A */ 01174 size_t ncol, /* # of columns of A */ 01175 size_t nzmax, /* max # of nonzeros of A */ 01176 int sorted, /* TRUE if columns of A sorted, FALSE otherwise */ 01177 int packed, /* TRUE if A will be packed, FALSE otherwise */ 01178 int stype, /* stype of A */ 01179 int xtype, /* CHOLMOD_PATTERN, _REAL, _COMPLEX, or _ZOMPLEX */ 01180 /* --------------- */ 01181 cholmod_common *Common 01182 ) ; 01183 01184 cholmod_sparse *cholmod_l_allocate_sparse (size_t, size_t, size_t, int, int, 01185 int, int, cholmod_common *) ; 01186 01187 /* -------------------------------------------------------------------------- */ 01188 /* cholmod_free_sparse: free a sparse matrix */ 01189 /* -------------------------------------------------------------------------- */ 01190 01191 int cholmod_free_sparse 01192 ( 01193 /* ---- in/out --- */ 01194 cholmod_sparse **A, /* matrix to deallocate, NULL on output */ 01195 /* --------------- */ 01196 cholmod_common *Common 01197 ) ; 01198 01199 int cholmod_l_free_sparse (cholmod_sparse **, cholmod_common *) ; 01200 01201 /* -------------------------------------------------------------------------- */ 01202 /* cholmod_reallocate_sparse: change the size (# entries) of sparse matrix */ 01203 /* -------------------------------------------------------------------------- */ 01204 01205 int cholmod_reallocate_sparse 01206 ( 01207 /* ---- input ---- */ 01208 size_t nznew, /* new # of entries in A */ 01209 /* ---- in/out --- */ 01210 cholmod_sparse *A, /* matrix to reallocate */ 01211 /* --------------- */ 01212 cholmod_common *Common 01213 ) ; 01214 01215 int cholmod_l_reallocate_sparse ( size_t, cholmod_sparse *, cholmod_common *) ; 01216 01217 /* -------------------------------------------------------------------------- */ 01218 /* cholmod_nnz: return number of nonzeros in a sparse matrix */ 01219 /* -------------------------------------------------------------------------- */ 01220 01221 UF_long cholmod_nnz 01222 ( 01223 /* ---- input ---- */ 01224 cholmod_sparse *A, 01225 /* --------------- */ 01226 cholmod_common *Common 01227 ) ; 01228 01229 UF_long cholmod_l_nnz (cholmod_sparse *, cholmod_common *) ; 01230 01231 /* -------------------------------------------------------------------------- */ 01232 /* cholmod_speye: sparse identity matrix */ 01233 /* -------------------------------------------------------------------------- */ 01234 01235 cholmod_sparse *cholmod_speye 01236 ( 01237 /* ---- input ---- */ 01238 size_t nrow, /* # of rows of A */ 01239 size_t ncol, /* # of columns of A */ 01240 int xtype, /* CHOLMOD_PATTERN, _REAL, _COMPLEX, or _ZOMPLEX */ 01241 /* --------------- */ 01242 cholmod_common *Common 01243 ) ; 01244 01245 cholmod_sparse *cholmod_l_speye (size_t, size_t, int, cholmod_common *) ; 01246 01247 /* -------------------------------------------------------------------------- */ 01248 /* cholmod_spzeros: sparse zero matrix */ 01249 /* -------------------------------------------------------------------------- */ 01250 01251 cholmod_sparse *cholmod_spzeros 01252 ( 01253 /* ---- input ---- */ 01254 size_t nrow, /* # of rows of A */ 01255 size_t ncol, /* # of columns of A */ 01256 size_t nzmax, /* max # of nonzeros of A */ 01257 int xtype, /* CHOLMOD_PATTERN, _REAL, _COMPLEX, or _ZOMPLEX */ 01258 /* --------------- */ 01259 cholmod_common *Common 01260 ) ; 01261 01262 cholmod_sparse *cholmod_l_spzeros (size_t, size_t, size_t, int, 01263 cholmod_common *) ; 01264 01265 /* -------------------------------------------------------------------------- */ 01266 /* cholmod_transpose: transpose a sparse matrix */ 01267 /* -------------------------------------------------------------------------- */ 01268 01269 /* Return A' or A.' The "values" parameter is 0, 1, or 2 to denote the pattern 01270 * transpose, the array transpose (A.'), and the complex conjugate transpose 01271 * (A'). 01272 */ 01273 01274 cholmod_sparse *cholmod_transpose 01275 ( 01276 /* ---- input ---- */ 01277 cholmod_sparse *A, /* matrix to transpose */ 01278 int values, /* 0: pattern, 1: array transpose, 2: conj. transpose */ 01279 /* --------------- */ 01280 cholmod_common *Common 01281 ) ; 01282 01283 cholmod_sparse *cholmod_l_transpose (cholmod_sparse *, int, cholmod_common *) ; 01284 01285 /* -------------------------------------------------------------------------- */ 01286 /* cholmod_transpose_unsym: transpose an unsymmetric sparse matrix */ 01287 /* -------------------------------------------------------------------------- */ 01288 01289 /* Compute F = A', A (:,f)', or A (p,f)', where A is unsymmetric and F is 01290 * already allocated. See cholmod_transpose for a simpler routine. */ 01291 01292 int cholmod_transpose_unsym 01293 ( 01294 /* ---- input ---- */ 01295 cholmod_sparse *A, /* matrix to transpose */ 01296 int values, /* 0: pattern, 1: array transpose, 2: conj. transpose */ 01297 int *Perm, /* size nrow, if present (can be NULL) */ 01298 int *fset, /* subset of 0:(A->ncol)-1 */ 01299 size_t fsize, /* size of fset */ 01300 /* ---- output --- */ 01301 cholmod_sparse *F, /* F = A', A(:,f)', or A(p,f)' */ 01302 /* --------------- */ 01303 cholmod_common *Common 01304 ) ; 01305 01306 int cholmod_l_transpose_unsym (cholmod_sparse *, int, UF_long *, UF_long *, 01307 size_t, cholmod_sparse *, cholmod_common *) ; 01308 01309 /* -------------------------------------------------------------------------- */ 01310 /* cholmod_transpose_sym: transpose a symmetric sparse matrix */ 01311 /* -------------------------------------------------------------------------- */ 01312 01313 /* Compute F = A' or A (p,p)', where A is symmetric and F is already allocated. 01314 * See cholmod_transpose for a simpler routine. */ 01315 01316 int cholmod_transpose_sym 01317 ( 01318 /* ---- input ---- */ 01319 cholmod_sparse *A, /* matrix to transpose */ 01320 int values, /* 0: pattern, 1: array transpose, 2: conj. transpose */ 01321 int *Perm, /* size nrow, if present (can be NULL) */ 01322 /* ---- output --- */ 01323 cholmod_sparse *F, /* F = A' or A(p,p)' */ 01324 /* --------------- */ 01325 cholmod_common *Common 01326 ) ; 01327 01328 int cholmod_l_transpose_sym (cholmod_sparse *, int, UF_long *, cholmod_sparse *, 01329 cholmod_common *) ; 01330 01331 /* -------------------------------------------------------------------------- */ 01332 /* cholmod_ptranspose: transpose a sparse matrix */ 01333 /* -------------------------------------------------------------------------- */ 01334 01335 /* Return A' or A(p,p)' if A is symmetric. Return A', A(:,f)', or A(p,f)' if 01336 * A is unsymmetric. */ 01337 01338 cholmod_sparse *cholmod_ptranspose 01339 ( 01340 /* ---- input ---- */ 01341 cholmod_sparse *A, /* matrix to transpose */ 01342 int values, /* 0: pattern, 1: array transpose, 2: conj. transpose */ 01343 int *Perm, /* if non-NULL, F = A(p,f) or A(p,p) */ 01344 int *fset, /* subset of 0:(A->ncol)-1 */ 01345 size_t fsize, /* size of fset */ 01346 /* --------------- */ 01347 cholmod_common *Common 01348 ) ; 01349 01350 cholmod_sparse *cholmod_l_ptranspose (cholmod_sparse *, int, UF_long *, 01351 UF_long *, size_t, cholmod_common *) ; 01352 01353 /* -------------------------------------------------------------------------- */ 01354 /* cholmod_sort: sort row indices in each column of sparse matrix */ 01355 /* -------------------------------------------------------------------------- */ 01356 01357 int cholmod_sort 01358 ( 01359 /* ---- in/out --- */ 01360 cholmod_sparse *A, /* matrix to sort */ 01361 /* --------------- */ 01362 cholmod_common *Common 01363 ) ; 01364 01365 int cholmod_l_sort (cholmod_sparse *, cholmod_common *) ; 01366 01367 /* -------------------------------------------------------------------------- */ 01368 /* cholmod_band: C = tril (triu (A,k1), k2) */ 01369 /* -------------------------------------------------------------------------- */ 01370 01371 cholmod_sparse *cholmod_band 01372 ( 01373 /* ---- input ---- */ 01374 cholmod_sparse *A, /* matrix to extract band matrix from */ 01375 UF_long k1, /* ignore entries below the k1-st diagonal */ 01376 UF_long k2, /* ignore entries above the k2-nd diagonal */ 01377 int mode, /* >0: numerical, 0: pattern, <0: pattern (no diag) */ 01378 /* --------------- */ 01379 cholmod_common *Common 01380 ) ; 01381 01382 cholmod_sparse *cholmod_l_band (cholmod_sparse *, UF_long, UF_long, int, 01383 cholmod_common *) ; 01384 01385 /* -------------------------------------------------------------------------- */ 01386 /* cholmod_band_inplace: A = tril (triu (A,k1), k2) */ 01387 /* -------------------------------------------------------------------------- */ 01388 01389 int cholmod_band_inplace 01390 ( 01391 /* ---- input ---- */ 01392 UF_long k1, /* ignore entries below the k1-st diagonal */ 01393 UF_long k2, /* ignore entries above the k2-nd diagonal */ 01394 int mode, /* >0: numerical, 0: pattern, <0: pattern (no diag) */ 01395 /* ---- in/out --- */ 01396 cholmod_sparse *A, /* matrix from which entries not in band are removed */ 01397 /* --------------- */ 01398 cholmod_common *Common 01399 ) ; 01400 01401 int cholmod_l_band_inplace (UF_long, UF_long, int, cholmod_sparse *, 01402 cholmod_common *) ; 01403 01404 /* -------------------------------------------------------------------------- */ 01405 /* cholmod_aat: C = A*A' or A(:,f)*A(:,f)' */ 01406 /* -------------------------------------------------------------------------- */ 01407 01408 cholmod_sparse *cholmod_aat 01409 ( 01410 /* ---- input ---- */ 01411 cholmod_sparse *A, /* input matrix; C=A*A' is constructed */ 01412 int *fset, /* subset of 0:(A->ncol)-1 */ 01413 size_t fsize, /* size of fset */ 01414 int mode, /* >0: numerical, 0: pattern, <0: pattern (no diag), 01415 * -2: pattern only, no diagonal, add 50%+n extra 01416 * space to C */ 01417 /* --------------- */ 01418 cholmod_common *Common 01419 ) ; 01420 01421 cholmod_sparse *cholmod_l_aat (cholmod_sparse *, UF_long *, size_t, int, 01422 cholmod_common *) ; 01423 01424 /* -------------------------------------------------------------------------- */ 01425 /* cholmod_copy_sparse: C = A, create an exact copy of a sparse matrix */ 01426 /* -------------------------------------------------------------------------- */ 01427 01428 cholmod_sparse *cholmod_copy_sparse 01429 ( 01430 /* ---- input ---- */ 01431 cholmod_sparse *A, /* matrix to copy */ 01432 /* --------------- */ 01433 cholmod_common *Common 01434 ) ; 01435 01436 cholmod_sparse *cholmod_l_copy_sparse (cholmod_sparse *, cholmod_common *) ; 01437 01438 /* -------------------------------------------------------------------------- */ 01439 /* cholmod_copy: C = A, with possible change of stype */ 01440 /* -------------------------------------------------------------------------- */ 01441 01442 cholmod_sparse *cholmod_copy 01443 ( 01444 /* ---- input ---- */ 01445 cholmod_sparse *A, /* matrix to copy */ 01446 int stype, /* requested stype of C */ 01447 int mode, /* >0: numerical, 0: pattern, <0: pattern (no diag) */ 01448 /* --------------- */ 01449 cholmod_common *Common 01450 ) ; 01451 01452 cholmod_sparse *cholmod_l_copy (cholmod_sparse *, int, int, cholmod_common *) ; 01453 01454 /* -------------------------------------------------------------------------- */ 01455 /* cholmod_add: C = alpha*A + beta*B */ 01456 /* -------------------------------------------------------------------------- */ 01457 01458 cholmod_sparse *cholmod_add 01459 ( 01460 /* ---- input ---- */ 01461 cholmod_sparse *A, /* matrix to add */ 01462 cholmod_sparse *B, /* matrix to add */ 01463 double alpha [2], /* scale factor for A */ 01464 double beta [2], /* scale factor for B */ 01465 int values, /* if TRUE compute the numerical values of C */ 01466 int sorted, /* if TRUE, sort columns of C */ 01467 /* --------------- */ 01468 cholmod_common *Common 01469 ) ; 01470 01471 cholmod_sparse *cholmod_l_add (cholmod_sparse *, cholmod_sparse *, double *, 01472 double *, int, int, cholmod_common *) ; 01473 01474 /* -------------------------------------------------------------------------- */ 01475 /* cholmod_sparse_xtype: change the xtype of a sparse matrix */ 01476 /* -------------------------------------------------------------------------- */ 01477 01478 int cholmod_sparse_xtype 01479 ( 01480 /* ---- input ---- */ 01481 int to_xtype, /* requested xtype (pattern, real, complex, zomplex) */ 01482 /* ---- in/out --- */ 01483 cholmod_sparse *A, /* sparse matrix to change */ 01484 /* --------------- */ 01485 cholmod_common *Common 01486 ) ; 01487 01488 int cholmod_l_sparse_xtype (int, cholmod_sparse *, cholmod_common *) ; 01489 01490 01491 /* ========================================================================== */ 01492 /* === Core/cholmod_factor ================================================== */ 01493 /* ========================================================================== */ 01494 01495 /* A symbolic and numeric factorization, either simplicial or supernodal. 01496 * In all cases, the row indices in the columns of L are kept sorted. */ 01497 01498 typedef struct cholmod_factor_struct 01499 { 01500 /* ---------------------------------------------------------------------- */ 01501 /* for both simplicial and supernodal factorizations */ 01502 /* ---------------------------------------------------------------------- */ 01503 01504 size_t n ; /* L is n-by-n */ 01505 01506 size_t minor ; /* If the factorization failed, L->minor is the column 01507 * at which it failed (in the range 0 to n-1). A value 01508 * of n means the factorization was successful or 01509 * the matrix has not yet been factorized. */ 01510 01511 /* ---------------------------------------------------------------------- */ 01512 /* symbolic ordering and analysis */ 01513 /* ---------------------------------------------------------------------- */ 01514 01515 void *Perm ; /* size n, permutation used */ 01516 void *ColCount ; /* size n, column counts for simplicial L */ 01517 01518 /* ---------------------------------------------------------------------- */ 01519 /* simplicial factorization */ 01520 /* ---------------------------------------------------------------------- */ 01521 01522 size_t nzmax ; /* size of i and x */ 01523 01524 void *p ; /* p [0..ncol], the column pointers */ 01525 void *i ; /* i [0..nzmax-1], the row indices */ 01526 void *x ; /* x [0..nzmax-1], the numerical values */ 01527 void *z ; 01528 void *nz ; /* nz [0..ncol-1], the # of nonzeros in each column. 01529 * i [p [j] ... p [j]+nz[j]-1] contains the row indices, 01530 * and the numerical values are in the same locatins 01531 * in x. The value of i [p [k]] is always k. */ 01532 01533 void *next ; /* size ncol+2. next [j] is the next column in i/x */ 01534 void *prev ; /* size ncol+2. prev [j] is the prior column in i/x. 01535 * head of the list is ncol+1, and the tail is ncol. */ 01536 01537 /* ---------------------------------------------------------------------- */ 01538 /* supernodal factorization */ 01539 /* ---------------------------------------------------------------------- */ 01540 01541 /* Note that L->x is shared with the simplicial data structure. L->x has 01542 * size L->nzmax for a simplicial factor, and size L->xsize for a supernodal 01543 * factor. */ 01544 01545 size_t nsuper ; /* number of supernodes */ 01546 size_t ssize ; /* size of s, integer part of supernodes */ 01547 size_t xsize ; /* size of x, real part of supernodes */ 01548 size_t maxcsize ; /* size of largest update matrix */ 01549 size_t maxesize ; /* max # of rows in supernodes, excl. triangular part */ 01550 01551 void *super ; /* size nsuper+1, first col in each supernode */ 01552 void *pi ; /* size nsuper+1, pointers to integer patterns */ 01553 void *px ; /* size nsuper+1, pointers to real parts */ 01554 void *s ; /* size ssize, integer part of supernodes */ 01555 01556 /* ---------------------------------------------------------------------- */ 01557 /* factorization type */ 01558 /* ---------------------------------------------------------------------- */ 01559 01560 int ordering ; /* ordering method used */ 01561 01562 int is_ll ; /* TRUE if LL', FALSE if LDL' */ 01563 int is_super ; /* TRUE if supernodal, FALSE if simplicial */ 01564 int is_monotonic ; /* TRUE if columns of L appear in order 0..n-1. 01565 * Only applicable to simplicial numeric types. */ 01566 01567 /* There are 8 types of factor objects that cholmod_factor can represent 01568 * (only 6 are used): 01569 * 01570 * Numeric types (xtype is not CHOLMOD_PATTERN) 01571 * -------------------------------------------- 01572 * 01573 * simplicial LDL': (is_ll FALSE, is_super FALSE). Stored in compressed 01574 * column form, using the simplicial components above (nzmax, p, i, 01575 * x, z, nz, next, and prev). The unit diagonal of L is not stored, 01576 * and D is stored in its place. There are no supernodes. 01577 * 01578 * simplicial LL': (is_ll TRUE, is_super FALSE). Uses the same storage 01579 * scheme as the simplicial LDL', except that D does not appear. 01580 * The first entry of each column of L is the diagonal entry of 01581 * that column of L. 01582 * 01583 * supernodal LDL': (is_ll FALSE, is_super TRUE). Not used. 01584 * FUTURE WORK: add support for supernodal LDL' 01585 * 01586 * supernodal LL': (is_ll TRUE, is_super TRUE). A supernodal factor, 01587 * using the supernodal components described above (nsuper, ssize, 01588 * xsize, maxcsize, maxesize, super, pi, px, s, x, and z). 01589 * 01590 * 01591 * Symbolic types (xtype is CHOLMOD_PATTERN) 01592 * ----------------------------------------- 01593 * 01594 * simplicial LDL': (is_ll FALSE, is_super FALSE). Nothing is present 01595 * except Perm and ColCount. 01596 * 01597 * simplicial LL': (is_ll TRUE, is_super FALSE). Identical to the 01598 * simplicial LDL', except for the is_ll flag. 01599 * 01600 * supernodal LDL': (is_ll FALSE, is_super TRUE). Not used. 01601 * FUTURE WORK: add support for supernodal LDL' 01602 * 01603 * supernodal LL': (is_ll TRUE, is_super TRUE). A supernodal symbolic 01604 * factorization. The simplicial symbolic information is present 01605 * (Perm and ColCount), as is all of the supernodal factorization 01606 * except for the numerical values (x and z). 01607 */ 01608 01609 int itype ; /* The integer arrays are Perm, ColCount, p, i, nz, 01610 * next, prev, super, pi, px, and s. If itype is 01611 * CHOLMOD_INT, all of these are int arrays. 01612 * CHOLMOD_INTLONG: p, pi, px are UF_long, others int. 01613 * CHOLMOD_LONG: all integer arrays are UF_long. */ 01614 int xtype ; /* pattern, real, complex, or zomplex */ 01615 int dtype ; /* x and z double or float */ 01616 01617 } cholmod_factor ; 01618 01619 01620 /* -------------------------------------------------------------------------- */ 01621 /* cholmod_allocate_factor: allocate a factor (symbolic LL' or LDL') */ 01622 /* -------------------------------------------------------------------------- */ 01623 01624 cholmod_factor *cholmod_allocate_factor 01625 ( 01626 /* ---- input ---- */ 01627 size_t n, /* L is n-by-n */ 01628 /* --------------- */ 01629 cholmod_common *Common 01630 ) ; 01631 01632 cholmod_factor *cholmod_l_allocate_factor (size_t, cholmod_common *) ; 01633 01634 /* -------------------------------------------------------------------------- */ 01635 /* cholmod_free_factor: free a factor */ 01636 /* -------------------------------------------------------------------------- */ 01637 01638 int cholmod_free_factor 01639 ( 01640 /* ---- in/out --- */ 01641 cholmod_factor **L, /* factor to free, NULL on output */ 01642 /* --------------- */ 01643 cholmod_common *Common 01644 ) ; 01645 01646 int cholmod_l_free_factor (cholmod_factor **, cholmod_common *) ; 01647 01648 /* -------------------------------------------------------------------------- */ 01649 /* cholmod_reallocate_factor: change the # entries in a factor */ 01650 /* -------------------------------------------------------------------------- */ 01651 01652 int cholmod_reallocate_factor 01653 ( 01654 /* ---- input ---- */ 01655 size_t nznew, /* new # of entries in L */ 01656 /* ---- in/out --- */ 01657 cholmod_factor *L, /* factor to modify */ 01658 /* --------------- */ 01659 cholmod_common *Common 01660 ) ; 01661 01662 int cholmod_l_reallocate_factor (size_t, cholmod_factor *, cholmod_common *) ; 01663 01664 /* -------------------------------------------------------------------------- */ 01665 /* cholmod_change_factor: change the type of factor (e.g., LDL' to LL') */ 01666 /* -------------------------------------------------------------------------- */ 01667 01668 int cholmod_change_factor 01669 ( 01670 /* ---- input ---- */ 01671 int to_xtype, /* to CHOLMOD_PATTERN, _REAL, _COMPLEX, _ZOMPLEX */ 01672 int to_ll, /* TRUE: convert to LL', FALSE: LDL' */ 01673 int to_super, /* TRUE: convert to supernodal, FALSE: simplicial */ 01674 int to_packed, /* TRUE: pack simplicial columns, FALSE: do not pack */ 01675 int to_monotonic, /* TRUE: put simplicial columns in order, FALSE: not */ 01676 /* ---- in/out --- */ 01677 cholmod_factor *L, /* factor to modify */ 01678 /* --------------- */ 01679 cholmod_common *Common 01680 ) ; 01681 01682 int cholmod_l_change_factor ( int, int, int, int, int, cholmod_factor *, 01683 cholmod_common *) ; 01684 01685 /* -------------------------------------------------------------------------- */ 01686 /* cholmod_pack_factor: pack the columns of a factor */ 01687 /* -------------------------------------------------------------------------- */ 01688 01689 /* Pack the columns of a simplicial factor. Unlike cholmod_change_factor, 01690 * it can pack the columns of a factor even if they are not stored in their 01691 * natural order (non-monotonic). */ 01692 01693 int cholmod_pack_factor 01694 ( 01695 /* ---- in/out --- */ 01696 cholmod_factor *L, /* factor to modify */ 01697 /* --------------- */ 01698 cholmod_common *Common 01699 ) ; 01700 01701 int cholmod_l_pack_factor (cholmod_factor *, cholmod_common *) ; 01702 01703 /* -------------------------------------------------------------------------- */ 01704 /* cholmod_reallocate_column: resize a single column of a factor */ 01705 /* -------------------------------------------------------------------------- */ 01706 01707 int cholmod_reallocate_column 01708 ( 01709 /* ---- input ---- */ 01710 size_t j, /* the column to reallocate */ 01711 size_t need, /* required size of column j */ 01712 /* ---- in/out --- */ 01713 cholmod_factor *L, /* factor to modify */ 01714 /* --------------- */ 01715 cholmod_common *Common 01716 ) ; 01717 01718 int cholmod_l_reallocate_column (size_t, size_t, cholmod_factor *, 01719 cholmod_common *) ; 01720 01721 /* -------------------------------------------------------------------------- */ 01722 /* cholmod_factor_to_sparse: create a sparse matrix copy of a factor */ 01723 /* -------------------------------------------------------------------------- */ 01724 01725 /* Only operates on numeric factors, not symbolic ones */ 01726 01727 cholmod_sparse *cholmod_factor_to_sparse 01728 ( 01729 /* ---- in/out --- */ 01730 cholmod_factor *L, /* factor to copy, converted to symbolic on output */ 01731 /* --------------- */ 01732 cholmod_common *Common 01733 ) ; 01734 01735 cholmod_sparse *cholmod_l_factor_to_sparse (cholmod_factor *, 01736 cholmod_common *) ; 01737 01738 /* -------------------------------------------------------------------------- */ 01739 /* cholmod_copy_factor: create a copy of a factor */ 01740 /* -------------------------------------------------------------------------- */ 01741 01742 cholmod_factor *cholmod_copy_factor 01743 ( 01744 /* ---- input ---- */ 01745 cholmod_factor *L, /* factor to copy */ 01746 /* --------------- */ 01747 cholmod_common *Common 01748 ) ; 01749 01750 cholmod_factor *cholmod_l_copy_factor (cholmod_factor *, cholmod_common *) ; 01751 01752 /* -------------------------------------------------------------------------- */ 01753 /* cholmod_factor_xtype: change the xtype of a factor */ 01754 /* -------------------------------------------------------------------------- */ 01755 01756 int cholmod_factor_xtype 01757 ( 01758 /* ---- input ---- */ 01759 int to_xtype, /* requested xtype (real, complex, or zomplex) */ 01760 /* ---- in/out --- */ 01761 cholmod_factor *L, /* factor to change */ 01762 /* --------------- */ 01763 cholmod_common *Common 01764 ) ; 01765 01766 int cholmod_l_factor_xtype (int, cholmod_factor *, cholmod_common *) ; 01767 01768 01769 /* ========================================================================== */ 01770 /* === Core/cholmod_dense =================================================== */ 01771 /* ========================================================================== */ 01772 01773 /* A dense matrix in column-oriented form. It has no itype since it contains 01774 * no integers. Entry in row i and column j is located in x [i+j*d]. 01775 */ 01776 01777 typedef struct cholmod_dense_struct 01778 { 01779 size_t nrow ; /* the matrix is nrow-by-ncol */ 01780 size_t ncol ; 01781 size_t nzmax ; /* maximum number of entries in the matrix */ 01782 size_t d ; /* leading dimension (d >= nrow must hold) */ 01783 void *x ; /* size nzmax or 2*nzmax, if present */ 01784 void *z ; /* size nzmax, if present */ 01785 int xtype ; /* pattern, real, complex, or zomplex */ 01786 int dtype ; /* x and z double or float */ 01787 01788 } cholmod_dense ; 01789 01790 /* -------------------------------------------------------------------------- */ 01791 /* cholmod_allocate_dense: allocate a dense matrix (contents uninitialized) */ 01792 /* -------------------------------------------------------------------------- */ 01793 01794 cholmod_dense *cholmod_allocate_dense 01795 ( 01796 /* ---- input ---- */ 01797 size_t nrow, /* # of rows of matrix */ 01798 size_t ncol, /* # of columns of matrix */ 01799 size_t d, /* leading dimension */ 01800 int xtype, /* CHOLMOD_REAL, _COMPLEX, or _ZOMPLEX */ 01801 /* --------------- */ 01802 cholmod_common *Common 01803 ) ; 01804 01805 cholmod_dense *cholmod_l_allocate_dense (size_t, size_t, size_t, int, 01806 cholmod_common *) ; 01807 01808 /* -------------------------------------------------------------------------- */ 01809 /* cholmod_zeros: allocate a dense matrix and set it to zero */ 01810 /* -------------------------------------------------------------------------- */ 01811 01812 cholmod_dense *cholmod_zeros 01813 ( 01814 /* ---- input ---- */ 01815 size_t nrow, /* # of rows of matrix */ 01816 size_t ncol, /* # of columns of matrix */ 01817 int xtype, /* CHOLMOD_REAL, _COMPLEX, or _ZOMPLEX */ 01818 /* --------------- */ 01819 cholmod_common *Common 01820 ) ; 01821 01822 cholmod_dense *cholmod_l_zeros (size_t, size_t, int, cholmod_common *) ; 01823 01824 /* -------------------------------------------------------------------------- */ 01825 /* cholmod_ones: allocate a dense matrix and set it to all ones */ 01826 /* -------------------------------------------------------------------------- */ 01827 01828 cholmod_dense *cholmod_ones 01829 ( 01830 /* ---- input ---- */ 01831 size_t nrow, /* # of rows of matrix */ 01832 size_t ncol, /* # of columns of matrix */ 01833 int xtype, /* CHOLMOD_REAL, _COMPLEX, or _ZOMPLEX */ 01834 /* --------------- */ 01835 cholmod_common *Common 01836 ) ; 01837 01838 cholmod_dense *cholmod_l_ones (size_t, size_t, int, cholmod_common *) ; 01839 01840 /* -------------------------------------------------------------------------- */ 01841 /* cholmod_eye: allocate a dense matrix and set it to the identity matrix */ 01842 /* -------------------------------------------------------------------------- */ 01843 01844 cholmod_dense *cholmod_eye 01845 ( 01846 /* ---- input ---- */ 01847 size_t nrow, /* # of rows of matrix */ 01848 size_t ncol, /* # of columns of matrix */ 01849 int xtype, /* CHOLMOD_REAL, _COMPLEX, or _ZOMPLEX */ 01850 /* --------------- */ 01851 cholmod_common *Common 01852 ) ; 01853 01854 cholmod_dense *cholmod_l_eye (size_t, size_t, int, cholmod_common *) ; 01855 01856 /* -------------------------------------------------------------------------- */ 01857 /* cholmod_free_dense: free a dense matrix */ 01858 /* -------------------------------------------------------------------------- */ 01859 01860 int cholmod_free_dense 01861 ( 01862 /* ---- in/out --- */ 01863 cholmod_dense **X, /* dense matrix to deallocate, NULL on output */ 01864 /* --------------- */ 01865 cholmod_common *Common 01866 ) ; 01867 01868 int cholmod_l_free_dense (cholmod_dense **, cholmod_common *) ; 01869 01870 /* -------------------------------------------------------------------------- */ 01871 /* cholmod_sparse_to_dense: create a dense matrix copy of a sparse matrix */ 01872 /* -------------------------------------------------------------------------- */ 01873 01874 cholmod_dense *cholmod_sparse_to_dense 01875 ( 01876 /* ---- input ---- */ 01877 cholmod_sparse *A, /* matrix to copy */ 01878 /* --------------- */ 01879 cholmod_common *Common 01880 ) ; 01881 01882 cholmod_dense *cholmod_l_sparse_to_dense (cholmod_sparse *, 01883 cholmod_common *) ; 01884 01885 /* -------------------------------------------------------------------------- */ 01886 /* cholmod_dense_to_sparse: create a sparse matrix copy of a dense matrix */ 01887 /* -------------------------------------------------------------------------- */ 01888 01889 cholmod_sparse *cholmod_dense_to_sparse 01890 ( 01891 /* ---- input ---- */ 01892 cholmod_dense *X, /* matrix to copy */ 01893 int values, /* TRUE if values to be copied, FALSE otherwise */ 01894 /* --------------- */ 01895 cholmod_common *Common 01896 ) ; 01897 01898 cholmod_sparse *cholmod_l_dense_to_sparse (cholmod_dense *, int, 01899 cholmod_common *) ; 01900 01901 /* -------------------------------------------------------------------------- */ 01902 /* cholmod_copy_dense: create a copy of a dense matrix */ 01903 /* -------------------------------------------------------------------------- */ 01904 01905 cholmod_dense *cholmod_copy_dense 01906 ( 01907 /* ---- input ---- */ 01908 cholmod_dense *X, /* matrix to copy */ 01909 /* --------------- */ 01910 cholmod_common *Common 01911 ) ; 01912 01913 cholmod_dense *cholmod_l_copy_dense (cholmod_dense *, cholmod_common *) ; 01914 01915 /* -------------------------------------------------------------------------- */ 01916 /* cholmod_copy_dense2: copy a dense matrix (pre-allocated) */ 01917 /* -------------------------------------------------------------------------- */ 01918 01919 int cholmod_copy_dense2 01920 ( 01921 /* ---- input ---- */ 01922 cholmod_dense *X, /* matrix to copy */ 01923 /* ---- output --- */ 01924 cholmod_dense *Y, /* copy of matrix X */ 01925 /* --------------- */ 01926 cholmod_common *Common 01927 ) ; 01928 01929 int cholmod_l_copy_dense2 (cholmod_dense *, cholmod_dense *, cholmod_common *) ; 01930 01931 /* -------------------------------------------------------------------------- */ 01932 /* cholmod_dense_xtype: change the xtype of a dense matrix */ 01933 /* -------------------------------------------------------------------------- */ 01934 01935 int cholmod_dense_xtype 01936 ( 01937 /* ---- input ---- */ 01938 int to_xtype, /* requested xtype (real, complex,or zomplex) */ 01939 /* ---- in/out --- */ 01940 cholmod_dense *X, /* dense matrix to change */ 01941 /* --------------- */ 01942 cholmod_common *Common 01943 ) ; 01944 01945 int cholmod_l_dense_xtype (int, cholmod_dense *, cholmod_common *) ; 01946 01947 01948 /* ========================================================================== */ 01949 /* === Core/cholmod_triplet ================================================= */ 01950 /* ========================================================================== */ 01951 01952 /* A sparse matrix stored in triplet form. */ 01953 01954 typedef struct cholmod_triplet_struct 01955 { 01956 size_t nrow ; /* the matrix is nrow-by-ncol */ 01957 size_t ncol ; 01958 size_t nzmax ; /* maximum number of entries in the matrix */ 01959 size_t nnz ; /* number of nonzeros in the matrix */ 01960 01961 void *i ; /* i [0..nzmax-1], the row indices */ 01962 void *j ; /* j [0..nzmax-1], the column indices */ 01963 void *x ; /* size nzmax or 2*nzmax, if present */ 01964 void *z ; /* size nzmax, if present */ 01965 01966 int stype ; /* Describes what parts of the matrix are considered: 01967 * 01968 * 0: matrix is "unsymmetric": use both upper and lower triangular parts 01969 * (the matrix may actually be symmetric in pattern and value, but 01970 * both parts are explicitly stored and used). May be square or 01971 * rectangular. 01972 * >0: matrix is square and symmetric. Entries in the lower triangular 01973 * part are transposed and added to the upper triangular part when 01974 * the matrix is converted to cholmod_sparse form. 01975 * <0: matrix is square and symmetric. Entries in the upper triangular 01976 * part are transposed and added to the lower triangular part when 01977 * the matrix is converted to cholmod_sparse form. 01978 * 01979 * Note that stype>0 and stype<0 are different for cholmod_sparse and 01980 * cholmod_triplet. The reason is simple. You can permute a symmetric 01981 * triplet matrix by simply replacing a row and column index with their 01982 * new row and column indices, via an inverse permutation. Suppose 01983 * P = L->Perm is your permutation, and Pinv is an array of size n. 01984 * Suppose a symmetric matrix A is represent by a triplet matrix T, with 01985 * entries only in the upper triangular part. Then the following code: 01986 * 01987 * Ti = T->i ; 01988 * Tj = T->j ; 01989 * for (k = 0 ; k < n ; k++) Pinv [P [k]] = k ; 01990 * for (k = 0 ; k < nz ; k++) Ti [k] = Pinv [Ti [k]] ; 01991 * for (k = 0 ; k < nz ; k++) Tj [k] = Pinv [Tj [k]] ; 01992 * 01993 * creates the triplet form of C=P*A*P'. However, if T initially 01994 * contains just the upper triangular entries (T->stype = 1), after 01995 * permutation it has entries in both the upper and lower triangular 01996 * parts. These entries should be transposed when constructing the 01997 * cholmod_sparse form of A, which is what cholmod_triplet_to_sparse 01998 * does. Thus: 01999 * 02000 * C = cholmod_triplet_to_sparse (T, 0, &Common) ; 02001 * 02002 * will return the matrix C = P*A*P'. 02003 * 02004 * Since the triplet matrix T is so simple to generate, it's quite easy 02005 * to remove entries that you do not want, prior to converting T to the 02006 * cholmod_sparse form. So if you include these entries in T, CHOLMOD 02007 * assumes that there must be a reason (such as the one above). Thus, 02008 * no entry in a triplet matrix is ever ignored. 02009 */ 02010 02011 int itype ; /* CHOLMOD_LONG: i and j are UF_long. Otherwise int. */ 02012 int xtype ; /* pattern, real, complex, or zomplex */ 02013 int dtype ; /* x and z are double or float */ 02014 02015 } cholmod_triplet ; 02016 02017 /* -------------------------------------------------------------------------- */ 02018 /* cholmod_allocate_triplet: allocate a triplet matrix */ 02019 /* -------------------------------------------------------------------------- */ 02020 02021 cholmod_triplet *cholmod_allocate_triplet 02022 ( 02023 /* ---- input ---- */ 02024 size_t nrow, /* # of rows of T */ 02025 size_t ncol, /* # of columns of T */ 02026 size_t nzmax, /* max # of nonzeros of T */ 02027 int stype, /* stype of T */ 02028 int xtype, /* CHOLMOD_PATTERN, _REAL, _COMPLEX, or _ZOMPLEX */ 02029 /* --------------- */ 02030 cholmod_common *Common 02031 ) ; 02032 02033 cholmod_triplet *cholmod_l_allocate_triplet (size_t, size_t, size_t, int, int, 02034 cholmod_common *) ; 02035 02036 /* -------------------------------------------------------------------------- */ 02037 /* cholmod_free_triplet: free a triplet matrix */ 02038 /* -------------------------------------------------------------------------- */ 02039 02040 int cholmod_free_triplet 02041 ( 02042 /* ---- in/out --- */ 02043 cholmod_triplet **T, /* triplet matrix to deallocate, NULL on output */ 02044 /* --------------- */ 02045 cholmod_common *Common 02046 ) ; 02047 02048 int cholmod_l_free_triplet (cholmod_triplet **, cholmod_common *) ; 02049 02050 /* -------------------------------------------------------------------------- */ 02051 /* cholmod_reallocate_triplet: change the # of entries in a triplet matrix */ 02052 /* -------------------------------------------------------------------------- */ 02053 02054 int cholmod_reallocate_triplet 02055 ( 02056 /* ---- input ---- */ 02057 size_t nznew, /* new # of entries in T */ 02058 /* ---- in/out --- */ 02059 cholmod_triplet *T, /* triplet matrix to modify */ 02060 /* --------------- */ 02061 cholmod_common *Common 02062 ) ; 02063 02064 int cholmod_l_reallocate_triplet (size_t, cholmod_triplet *, cholmod_common *) ; 02065 02066 /* -------------------------------------------------------------------------- */ 02067 /* cholmod_sparse_to_triplet: create a triplet matrix copy of a sparse matrix*/ 02068 /* -------------------------------------------------------------------------- */ 02069 02070 cholmod_triplet *cholmod_sparse_to_triplet 02071 ( 02072 /* ---- input ---- */ 02073 cholmod_sparse *A, /* matrix to copy */ 02074 /* --------------- */ 02075 cholmod_common *Common 02076 ) ; 02077 02078 cholmod_triplet *cholmod_l_sparse_to_triplet (cholmod_sparse *, 02079 cholmod_common *) ; 02080 02081 /* -------------------------------------------------------------------------- */ 02082 /* cholmod_triplet_to_sparse: create a sparse matrix copy of a triplet matrix*/ 02083 /* -------------------------------------------------------------------------- */ 02084 02085 cholmod_sparse *cholmod_triplet_to_sparse 02086 ( 02087 /* ---- input ---- */ 02088 cholmod_triplet *T, /* matrix to copy */ 02089 size_t nzmax, /* allocate at least this much space in output matrix */ 02090 /* --------------- */ 02091 cholmod_common *Common 02092 ) ; 02093 02094 cholmod_sparse *cholmod_l_triplet_to_sparse (cholmod_triplet *, size_t, 02095 cholmod_common *) ; 02096 02097 /* -------------------------------------------------------------------------- */ 02098 /* cholmod_copy_triplet: create a copy of a triplet matrix */ 02099 /* -------------------------------------------------------------------------- */ 02100 02101 cholmod_triplet *cholmod_copy_triplet 02102 ( 02103 /* ---- input ---- */ 02104 cholmod_triplet *T, /* matrix to copy */ 02105 /* --------------- */ 02106 cholmod_common *Common 02107 ) ; 02108 02109 cholmod_triplet *cholmod_l_copy_triplet (cholmod_triplet *, cholmod_common *) ; 02110 02111 /* -------------------------------------------------------------------------- */ 02112 /* cholmod_triplet_xtype: change the xtype of a triplet matrix */ 02113 /* -------------------------------------------------------------------------- */ 02114 02115 int cholmod_triplet_xtype 02116 ( 02117 /* ---- input ---- */ 02118 int to_xtype, /* requested xtype (pattern, real, complex,or zomplex)*/ 02119 /* ---- in/out --- */ 02120 cholmod_triplet *T, /* triplet matrix to change */ 02121 /* --------------- */ 02122 cholmod_common *Common 02123 ) ; 02124 02125 int cholmod_l_triplet_xtype (int, cholmod_triplet *, cholmod_common *) ; 02126 02127 02128 /* ========================================================================== */ 02129 /* === Core/cholmod_memory ================================================== */ 02130 /* ========================================================================== */ 02131 02132 /* The user may make use of these, just like malloc and free. You can even 02133 * malloc an object and safely free it with cholmod_free, and visa versa 02134 * (except that the memory usage statistics will be corrupted). These routines 02135 * do differ from malloc and free. If cholmod_free is given a NULL pointer, 02136 * for example, it does nothing (unlike the ANSI free). cholmod_realloc does 02137 * not return NULL if given a non-NULL pointer and a nonzero size, even if it 02138 * fails (it returns the original pointer and sets an error code in 02139 * Common->status instead). 02140 * 02141 * CHOLMOD keeps track of the amount of memory it has allocated, and so the 02142 * cholmod_free routine also takes the size of the object being freed. This 02143 * is only used for statistics. If you, the user of CHOLMOD, pass the wrong 02144 * size, the only consequence is that the memory usage statistics will be 02145 * corrupted. 02146 */ 02147 02148 void *cholmod_malloc /* returns pointer to the newly malloc'd block */ 02149 ( 02150 /* ---- input ---- */ 02151 size_t n, /* number of items */ 02152 size_t size, /* size of each item */ 02153 /* --------------- */ 02154 cholmod_common *Common 02155 ) ; 02156 02157 void *cholmod_l_malloc (size_t, size_t, cholmod_common *) ; 02158 02159 void *cholmod_calloc /* returns pointer to the newly calloc'd block */ 02160 ( 02161 /* ---- input ---- */ 02162 size_t n, /* number of items */ 02163 size_t size, /* size of each item */ 02164 /* --------------- */ 02165 cholmod_common *Common 02166 ) ; 02167 02168 void *cholmod_l_calloc (size_t, size_t, cholmod_common *) ; 02169 02170 void *cholmod_free /* always returns NULL */ 02171 ( 02172 /* ---- input ---- */ 02173 size_t n, /* number of items */ 02174 size_t size, /* size of each item */ 02175 /* ---- in/out --- */ 02176 void *p, /* block of memory to free */ 02177 /* --------------- */ 02178 cholmod_common *Common 02179 ) ; 02180 02181 void *cholmod_l_free (size_t, size_t, void *, cholmod_common *) ; 02182 02183 void *cholmod_realloc /* returns pointer to reallocated block */ 02184 ( 02185 /* ---- input ---- */ 02186 size_t nnew, /* requested # of items in reallocated block */ 02187 size_t size, /* size of each item */ 02188 /* ---- in/out --- */ 02189 void *p, /* block of memory to realloc */ 02190 size_t *n, /* current size on input, nnew on output if successful*/ 02191 /* --------------- */ 02192 cholmod_common *Common 02193 ) ; 02194 02195 void *cholmod_l_realloc (size_t, size_t, void *, size_t *, cholmod_common *) ; 02196 02197 int cholmod_realloc_multiple 02198 ( 02199 /* ---- input ---- */ 02200 size_t nnew, /* requested # of items in reallocated blocks */ 02201 int nint, /* number of int/UF_long blocks */ 02202 int xtype, /* CHOLMOD_PATTERN, _REAL, _COMPLEX, or _ZOMPLEX */ 02203 /* ---- in/out --- */ 02204 void **I, /* int or UF_long block */ 02205 void **J, /* int or UF_long block */ 02206 void **X, /* complex, double, or float block */ 02207 void **Z, /* zomplex case only: double or float block */ 02208 size_t *n, /* current size of the I,J,X,Z blocks on input, 02209 * nnew on output if successful */ 02210 /* --------------- */ 02211 cholmod_common *Common 02212 ) ; 02213 02214 int cholmod_l_realloc_multiple (size_t, int, int, void **, void **, void **, 02215 void **, size_t *, cholmod_common *) ; 02216 02217 /* ========================================================================== */ 02218 /* === symmetry types ======================================================= */ 02219 /* ========================================================================== */ 02220 02221 #define CHOLMOD_MM_RECTANGULAR 1 02222 #define CHOLMOD_MM_UNSYMMETRIC 2 02223 #define CHOLMOD_MM_SYMMETRIC 3 02224 #define CHOLMOD_MM_HERMITIAN 4 02225 #define CHOLMOD_MM_SKEW_SYMMETRIC 5 02226 #define CHOLMOD_MM_SYMMETRIC_POSDIAG 6 02227 #define CHOLMOD_MM_HERMITIAN_POSDIAG 7 02228 02229 /* ========================================================================== */ 02230 /* === Numerical relop macros =============================================== */ 02231 /* ========================================================================== */ 02232 02233 /* These macros correctly handle the NaN case. 02234 * 02235 * CHOLMOD_IS_NAN(x): 02236 * True if x is NaN. False otherwise. The commonly-existing isnan(x) 02237 * function could be used, but it's not in Kernighan & Ritchie 2nd edition 02238 * (ANSI C89). It may appear in <math.h>, but I'm not certain about 02239 * portability. The expression x != x is true if and only if x is NaN, 02240 * according to the IEEE 754 floating-point standard. 02241 * 02242 * CHOLMOD_IS_ZERO(x): 02243 * True if x is zero. False if x is nonzero, NaN, or +/- Inf. 02244 * This is (x == 0) if the compiler is IEEE 754 compliant. 02245 * 02246 * CHOLMOD_IS_NONZERO(x): 02247 * True if x is nonzero, NaN, or +/- Inf. False if x zero. 02248 * This is (x != 0) if the compiler is IEEE 754 compliant. 02249 * 02250 * CHOLMOD_IS_LT_ZERO(x): 02251 * True if x is < zero or -Inf. False if x is >= 0, NaN, or +Inf. 02252 * This is (x < 0) if the compiler is IEEE 754 compliant. 02253 * 02254 * CHOLMOD_IS_GT_ZERO(x): 02255 * True if x is > zero or +Inf. False if x is <= 0, NaN, or -Inf. 02256 * This is (x > 0) if the compiler is IEEE 754 compliant. 02257 * 02258 * CHOLMOD_IS_LE_ZERO(x): 02259 * True if x is <= zero or -Inf. False if x is > 0, NaN, or +Inf. 02260 * This is (x <= 0) if the compiler is IEEE 754 compliant. 02261 */ 02262 02263 #ifdef CHOLMOD_WINDOWS 02264 02265 /* Yes, this is exceedingly ugly. Blame Microsoft, which hopelessly */ 02266 /* violates the IEEE 754 floating-point standard in a bizarre way. */ 02267 /* If you're using an IEEE 754-compliant compiler, then x != x is true */ 02268 /* iff x is NaN. For Microsoft, (x < x) is true iff x is NaN. */ 02269 /* So either way, this macro safely detects a NaN. */ 02270 #define CHOLMOD_IS_NAN(x) (((x) != (x)) || (((x) < (x)))) 02271 #define CHOLMOD_IS_ZERO(x) (((x) == 0.) && !CHOLMOD_IS_NAN(x)) 02272 #define CHOLMOD_IS_NONZERO(x) (((x) != 0.) || CHOLMOD_IS_NAN(x)) 02273 #define CHOLMOD_IS_LT_ZERO(x) (((x) < 0.) && !CHOLMOD_IS_NAN(x)) 02274 #define CHOLMOD_IS_GT_ZERO(x) (((x) > 0.) && !CHOLMOD_IS_NAN(x)) 02275 #define CHOLMOD_IS_LE_ZERO(x) (((x) <= 0.) && !CHOLMOD_IS_NAN(x)) 02276 02277 #else 02278 02279 /* These all work properly, according to the IEEE 754 standard ... except on */ 02280 /* a PC with windows. Works fine in Linux on the same PC... */ 02281 #define CHOLMOD_IS_NAN(x) ((x) != (x)) 02282 #define CHOLMOD_IS_ZERO(x) ((x) == 0.) 02283 #define CHOLMOD_IS_NONZERO(x) ((x) != 0.) 02284 #define CHOLMOD_IS_LT_ZERO(x) ((x) < 0.) 02285 #define CHOLMOD_IS_GT_ZERO(x) ((x) > 0.) 02286 #define CHOLMOD_IS_LE_ZERO(x) ((x) <= 0.) 02287 02288 #endif 02289 02290 #endif