Index of OperatorsopBernoulli
opBernoulli Bernoulli-ensemble operator. opBernoulli(M,N) creates an M-by-N Bernoulli ensemble, a matrix with iid +1/-1 entries. opBernoulli(M) creates a square M-by-M Bernoulli ensemble. opBernoulli(M,N,MODE) is the same as above, except that the parameter MODE controls the type of ensemble that is generated. The default is MODE=0 unless the overall memory requred exceeds 50 MBs. MODE = 0 (default): generates an explicit unnormalized matrix with random +1/-1 entries. The overall storage is O(M*N). MODE = 1: generates columns of the unnormalized matrix as the operator is applied. This allows for much larger ensembles since the matrix is implicit. The overall storage is O(M). MODE = 2: generates a scaled explicit matrix with unit-norm columns. MODE = 3: same as MODE=2, but the matrix is implicit (see MODE=1). Available operator properties: .mode gives the mode used to create the operator. opBinary
opBinary Binary (0/1) ensemble. opBinary(M,N) creates an M-by-N binary-ensemble operator. opBinary(M) creates a square M-by-M binary-ensemble. opGaussian(M,N,MODE) is the same as above, except that the parameter MODE controls the type of ensemble that is generated. The default is MODE=0 unless the overall memory requred exceeds 50 MBs. MODE = 0 (default): generates an explicit matrix with O(M*N) storage. MODE = 1: generates columns of the matrix as the operator is applied. This allows for much larger ensembles because the matrix is stored implicitly. The overall storage is O(M). opBlockDiag
opBlockDiag Operator-diagonal operator. B = opBlockDiag(OP1, OP2,...,OPN,OVERLAP) creates a compound block operator with the input operators OP1, OP2,... on the diagonal of B, e.g., B = DIAG([OP1 OP2 ... OPN]). When OVERLAP is a positive integer the blocks will be offset OVERLAP rows relative to the previous operator, when OVERLAP is negative the operators are offset by the absolute value of OVERLAP in columns. Note that choosing OVERLAP larger than the operator size may cause the matrix to become block antidiagonal. B = opBlockDiag(WEIGHT,OP1,...,OPN,OVERLAP) additionally weights each block by the elements of the vector WEIGHT. If only a single operator is given it is replicated as many times as there are weights. B = opBlockDiag(N,OP,OVERLAP) similar as above with WEIGHT equal to ones(N,1). This will cause operator OP to be repeated N times. opBlockOp
opBlockOp Blockwise application of operator on matrices. B = opBlockOp(M,N,OPIN,BR1,BC1,BR2,BC2) creates an operator that applies the given OPIN operator on two-dimensional data in a blockwise fashion. In the forward mode this means that the input vector is reshaped into an M-by-N matrix, which is then divided into blocks of size BR1-by-BC1. Next, we apply OPIN to each (vectorized) block and reshape the output to BR2-by-BC2 blocks. These blocks are gathered in a matrix which is vectorized to give the final output. In transpose mode, the input vector is reshaped into a matrix with M/BR1-by-N/BC1 blocks of size BR2-by-BC2, and the conjugate transpose of OPIN is applied to each block as described above to give BR1-by-BC1 blocks. These form an M-by-N matrix which is vectorized for output. When omitted, BR2 and BC2 are respectively set to BR1 and BC1 by default. opCTranspose
opCTranspose Conjugate transpose of an operator. opCTranspose(OP) returns the conjugate tranpose of OP. opChol
opCHOL Operator representing the Cholesky factorization of a
symmetric and definite matrix with optional iterative
refinement. Only the lower triangle of the input matrix
is referenced.
opChol(A) creates an operator for multiplication by the
inverse of the matrix A implicitly represented by its Cholesky
factorization. Optionally, iterative refinement is performed.
Note that A is an explicit matrix.
The following attributes may be changed by the user:
* nitref : the maximum number of iterative refinement steps (3)
* itref_tol : iterative refinement tolerance (1.0e-8)
* force_itref : force iterative refinement (false)
Dominique Orban <dominique.orban@gerad.ca>, 2014.
Copyright 2009, Ewout van den Berg and Michael P. Friedlander
See the file COPYING.txt for full copyright information.
Use the command 'spot.gpl' to locate this file. opClass
opClass Wrapper for classes. opClass(M,N,OBJ,CFLAG,LINFLAG) creates an M by N wrapper operator for the class instance OBJ. The only requirement on OBJ is that it implements the `mtimes' method. Optional arguments CFLAG and LINFLAG indicate whether the class implements a complex or real operator and whether it is linear or not. By default these fields are set to CFLAG=0, LINFLAG=1. opConj
opConj Take the elementwise conjugate of a complex operator. opConj(OP) is the elementwise complex conjugate of operator OP. Applying opConj to conjugate operators returns the original operator. opConvolve
opConvolve One and two dimensional convolution operator. opConvolve(M,N,KERNEL,OFFSET,MODE) creates an operator for one or two-dimensional convolution, depending on the size of the KERNEL, and the matrix or vector (MxN) the convolution is applied to. The convolution is one dimensional only if KERNEL is a column vector and N=1, or KERNEL is a row vector and M=1. The OFFSET parameter determines the center of the KERNEL and has a default value of [1,1]. When the OFFSET lies outside the size of the KERNEL, the KERNEL is embedded in a zero matrix/vector with appropriate center. For one-dimensional convolution, KERNEL may be a scalar. Specifying an offset that is not equal to one where the corresponding size of the kernel does equal one leads to the construction of a two-dimensional convolution operator. There are three types of MODE: MODE = 'regular' - convolve input with kernel; 'truncated' - convolve input with kernel, but keep only those MxN entries in the result that overlap with the input; 'cyclic' - do cyclic convolution of the input with a kernel that is wrapped around as many times as needed. The output of the convolution operator, like all other operators, is in vector form. opCurvelet
opCurvelet Two-dimensional curvelet operator. opCurvelet(M,N,NBSCALES,NBANGLES,TTYPE) creates a two-dimensional curvelet operator for M by N matrices. The curvelet transform is computed using the Curvelab code. The remaining three parameters are optional; NBSCALES gives the number of scales and is set to max(1,ceil(log2(min(M,N)) - 3)) by default, as suggested by Curvelab. NBANGLES gives the number of angles at the second coarsest level which must be a multiple of four with a minimum of 8. By default NBANGLES is set to 16. TTYPE determines the type of transformation and is set to 'WRAP' by default. opDCT
opDCT Discrete cosine transform (DCT). opDCT(M) creates a one-dimensional discrete cosine transform operator for vectors of length M. opDCT2
opDCT2 Two-dimensional discrete cosine transform (DCT). opDCT2(M,N) creates a two-dimensional discrete cosine transform operator for matrices of size M-by-N. Input and output of the matrices is done in vectorized form. When N is omitted it is set to M by default. opDFT
opDFT Fast Fourier transform (DFT). opDFT(M) create a unitary one-dimensional discrete Fourier transform (DFT) for vectors of length M. opDFT(M,CENTERED), with the CENTERED flag set to true, creates a unitary DFT that shifts the zero-frequency component to the center of the spectrum. opDFT2
opDFT2 Two-dimensional fast Fourier transform (DFT). opDFT2(M,N) creates a two-dimensional normalized Fourier transform operator for matrices of size M by N. Input and output of the matrices is done in vectorized form. opDFT2(M,N,CENTERED) just like opDFT2(M,N), but with components shifted to have to zero-frequency component in the center of the spectrum, if the CENTERED flag is set to true. opDiag
opDiag Diagonal operator. opDiag(D) creates an operator for multiplication by the diagonal matrix that has a vector D on its diagonal. opDictionary
opDictionary Dictionary of concatenated operators. D = opDictionary(OP1,OP2,...OPn) creates a dictionary operator consisting of the concatenation of all operators, i.e., D = [ OP1, OP2, ..., OPn ]. In general, it's best to use Matlab's horizonal concatenation operations instead of calling opDictionary. (The two are equivalent.) opDirac
opDirac Dirac basis. opDirac(N) creates the square N-by-N identity operator. Without any arguments an operator corresponding to the scalar 1 is created. opEmpty
opEmpty Operator equivalent to empty matrix. opEmpty(M,N) creates an operator corresponding to an empty M-by-N matrix. At least one of M and N must be zero. opExcise
opExcise Excise rows or columns of an operator. opExcise(OP,IDX,TYPE) excises the entries in the rows or columns given by IDX, depending on whether TYPE = 'rows', or 'cols'. opExtend
opExtend Symmetric extension operator. opExtend(P,Q,PEXT,QEXT) creates an extension operator that acts on a "vectorized" matrix and appends a mirror symmetric boundary to the right and bottom portions of matrix. The original matrix is of size PxQ. The operator generates a matrix that is PEXTxQEXT matrix. The adjoint of the operator creates a PxQ matrix whose entries adjacent to the extension border are twice the value of the entries of the original matrix. Q QEXT *************++++++ * * + * * + * * + P ************* + + + PEXT +++++++++++++++++++ Example 1. Extend a 2-by-3 matrix into a 4-by-6 matrix: A = [1 2 3; 4 5 6]; E = opExtend(2, 3, 4, 6); reshape(E*A(:),4,6) Example 2. Requires imaging toolbox: I = double(imread('cameraman.tif')); E = opExtend(256,256,650,1400); Isup = reshape(E*I(:),650,1400); figure; imshow(uint8(Isup)); opEye
opEye Identity operator. opEye(M) creates the M-by-M identity operator. opEye(M,N) creates the M-by-N identity operator. If N is omitted it is set to M by default. Without any arguments an operator corresponding to the scalar 1 is created. opEye([M N]) is the same as the above. opFactorization
opFactorization Operator representing an inverse by way of a factorization.
Useful to avoid multiple factorizations as in opInverse
and save time by performing only forward and backsolves. opFoG
opFoG Forms the product of two operators. opFoG(OP1,OP2) creates an operator that successively applies each of the operators OP1, OP2 on a given input vector. In non-adjoint mode this is done in reverse order. The inputs must be either Spot operators or explicit Matlab matrices (including scalars). opFunction
opFunction Wrapper for functions. opFunction(M,N,FUN) creates a wrapper for function FUN, which corresponds to an M-by-N operator. The FUN parameter can be one of two types: 1) A handle to a function of the form FUN(X,MODE), where the operator is applied to X when MODE = 1, and the transpose is applied when MODE = 2; 2) A cell array of two function handles: {FUN,FUN_TRANSPOSE}, each of which requires only one parameter, X. opFunction(M,N,FUN,CFLAG,LINFLAG) additionally allows the arguments CFLAG and LINFLAG to indicate whether the function implements a complex or real operator and whether it is linear or not. The default values are CFLAG=0, LINFLAG=1. opGaussian
opGaussian Gaussian ensemble. opGaussian(M,N) creates an M-by-N Gaussian-ensemble operator. opGaussian(M) creates a square M-by-M Gaussian-ensemble. opGaussian(M,N,MODE) is the same as above, except that the parameter MODE controls the type of ensemble that is generated. The default is MODE=0 unless the overall memory requred exceeds 50 MBs. MODE = 0 (default): generates an explicit unnormalized matrix from the Normal distribution. The overall storage is O(M*N). MODE = 1: generates columns of the unnormalized matrix as the operator is applied. This allows for much larger ensembles because the matrix is implicit. The overall storage is O(M). MODE = 2: generates a explicit matrix with unit-norm columns. MODE = 3: same as MODE=2, but the matrix is implicit (see MODE=1). MODE = 4: generates an explicit matrix with orthonormal rows. This mode requires M <= N. Available operator properties: .mode is the mode used to create the operator. .seed is the seed used to initialize the RNG. opHaar
opHaar Haar wavelet. opHaar(N) creates a Haar Wavelet operator for 1-D signals of length N using 5 levels. N must be a power of 2. opHaar(N,LEVELS) optionally allows the number of LEVELS to be specified. opHaar(N,LEVELS,REDUNDANT) optionally specifies the boolean field REDUNDANT (default false). (See opWavelet for a description of this option.) opHaar2
opHaar2 2-D Haar Wavelet. opHaar2(M,N) creates a Haar Wavelet operator for 2-D signals of size M-by-N using 5 levels. M*N must be a power of 2. opHaar2(M,N,LEVELS) optionally allows the number of LEVELS to be specified. opHaar2(M,N,LEVELS,REDUNDANT) optionally specifies the boolean field REDUNDANT (default false). (See opWavelet for a description of this option.) opHadamard
opHadamard Hadamard matrix. opHadamard(N) creates a Hadamard operator for vectors of length N, where N is a power of two. Multiplication is done using a fast routine. opHadamard(N,NORMALIZED) is the same as above, except that the columns are scaled to unit two-norm. By default, the NORMALIZED flag is set to FALSE. opHeaviside
opHeaviside Heaviside operator. opHeaviside(N,NORMALIZED) creates an operator for multiplication by an N by N Heaviside matrix. These matrices have ones below and on the diagonal and zeros elsewhere. NORMALIZED is a flag indicating whether the columns should be scaled to unit Euclidean norm. By default the columns are unnormalized. opHermitian
opHermitian Convert a numeric matrix stored as its lower triangle into a Spot operator. opHermitian(A,DESCRIPTION) creates a hermitian operator that performs matrix-vector multiplication with matrix A. Only the lower triangle of the input matrix need to be stored. The optional parameter DESCRIPTION can be used to override the default operator name when printed. opImag
opImag Complex imaginary part of operator. opImag(OP) is the complex imaginary part of operator OP. Note that the resulting operator is real. opInverse
opInverse (Pseudo) inverse of operator. Ainv = opInverse(A) creates the (pseudo) inverse of a square operator. The product Ainv*b is then equivalent to A\b. opKron
opKron Kronecker tensor product. opKron(OP1,OP2,...OPn) creates an operator that is the Kronecker tensor product of OP1, OP2, ..., OPn. opLDL
opLDL Operator representing the LDL factorization of a symmetric matrix with optional iterative refinement. Only the lower triangle of the input matrix is referenced. This is currently only available for real matrices. If the matrix is known to be definite, opChol will be more efficient. opLDL(A) creates an operator for multiplication by the inverse of the matrix A implicitly represented by its LDL factorization. Optionally, iterative refinement is performed. Note that A is an explicit matrix. The following attributes may be changed by the user: * nitref : the maximum number of iterative refinement steps (3) * itref_tol : iterative refinement tolerance (1.0e-8) * force_itref : force iterative refinement (false) Dominique Orban <dominique.orban@gerad.ca>, 2014. Copyright 2009, Ewout van den Berg and Michael P. Friedlander See the file COPYING.txt for full copyright information. Use the command 'spot.gpl' to locate this file. opLU
opLU Operator representing the LU factorization of a matrix with optional iterative refinement. If the matrix is known to be symmetric, opLDL will be more efficient. If in addition, the matrix is known to be definite, opChol will be more efficient. opLU(A) creates an operator for multiplication by the (pseudo-) inverse of the matrix A implicitly represented by its LU factorization. Optionally, iterative refinement is performed. Note that A is an explicit matrix. The following attributes may be changed by the user: * nitref : the maximum number of iterative refinement steps (3) * itref_tol : iterative refinement tolerance (1.0e-8) * force_itref : force iterative refinement (false) Dominique Orban <dominique.orban@gerad.ca>, 2014. Copyright 2009, Ewout van den Berg and Michael P. Friedlander See the file COPYING.txt for full copyright information. Use the command 'spot.gpl' to locate this file. opMask
opMask Selection mask. opMask(N,IDX) creates a diagonal N-by-N operator that has ones only on those locations indicated by IDX. opMask(IDX) is the same as opMask(numel(IDX),IDX) when parameter IDX is logical. opMatrix
opMatrix Convert a numeric matrix into a Spot operator. opMatrix(A,DESCRIPTION) creates an operator that performs matrix-vector multiplication with matrix A. The optional parameter DESCRIPTION can be used to override the default operator name when printed. opMinus
opMinus Difference of two operators. opMinus(OP1,OP2) creates a compound operator representing OP1-OP2. opOnes
opOnes Operator equivalent to ones function. opOnes(M,N) creates an operator corresponding to an M by N matrix of ones. If parameter N is omitted it is set to M. opOrthogonal
opOrthogonal Abstract class for orthogonal operators. opOrthogonal methods: opOrthogonal - constructor mldivide - solves Ax=b via x=A'b. opPInverse
opPInverse Pseudo inverse of operator. Apinv = opPInverse(A) creates the pseudo inverse of a M-by-N operator A. The product Apinv*b is then equivalent to A\b. opPermutation
opPermutation Permutation operator. P = opPermutation(p) creates a permutation operator from a permutation vector. The product P*b is then equivalent to b(p). Dominique Orban <dominique.orban@gerad.ca>, 2014. opPower
opPower Raise operator to integer power. opPower(OP,P) creates the operator OP^P for integer values of P. When P = 0, the identity matrix is returned. When P < 0 we reformulate the operator as inv(OP^|P|). opQR
opQR Operator representing the QR factorization of a matrix with optional iterative refinement. If the matrix opQR(A) creates an operator for multiplication by the (pseudo-) inverse of the matrix A implicitly represented by its QR factorization. Optionally, iterative refinement is performed. Note that A is an explicit matrix. The following attributes may be changed by the user: * nitref : the maximum number of iterative refinement steps (3) * itref_tol : iterative refinement tolerance (1.0e-8) * force_itref : force iterative refinement (false) Dominique Orban <dominique.orban@gerad.ca>, 2014. Copyright 2009, Ewout van den Berg and Michael P. Friedlander See the file COPYING.txt for full copyright information. Use the command 'spot.gpl' to locate this file. opReal
opReal Real part of operator. opReal(OP) is the real part of operator OP. opRestriction
opRestriction Restriction operator. opRestriction(N,IDX) creates a restriction operator that selects the entries listed in the index vector IDX from an input vector of length N. The adjoint of the operator creates a zero vector of length N and fills the entries given by IDX with the input data. Algebraically, opRestriction(N,IDX) is equivalent to a matrix of size length(IDX)-by-N, where row i has a single 1 in column IDX(i). opSparseBinary
opSparseBinary Random sparse binary matrix. opSparseBinary(M,N) creates an M-by-N sparse binary matrix with min(M,8) nonzeros in each column. opSparseBinary(M,N,D) is the same as above, except that each column has min(M,D) nonzero elements. This matrix was suggested by Radu Berinde and Piotr Indyk, "Sparse recovery using sparse random matrices", MIT CSAIL TR 2008-001, January 2008, http://hdl.handle.net/1721.1/40089 Note that opSparseBinary calls RANDPERM and thus changes the state of RAND. opStack
opStack Stack of vertically concatenated operators. opStack(WEIGHTS, OP1, OP2, ...) creates a stacked operator consisting of the vertical concatenation of all operators; [WEIGHT1*OP1 WEIGHT2*OP2 ... WEIGHTn*OPn] If the same weight is to be applied to each operator, set WEIGHTS to a scalar. When WEIGHTS is empty [], it is set to one. The WEIGHT parameter can be omitted as long as OP1 is not a vector of length (n-1); in which case there is no way to decide whether it is a weight vector or operator. opSubsAsgn
opSubsAsgn Redefine rectangular subset of operator.
opSubsAsign(A,ROWIDX,COLIDX,B), index = ':' is valid. Size of B
must match that of rowidx and colidx. THe indices can exceed the
size of A (even if idx is given as logical) but must not be negative? opSubsRef
opSubsRef Extract rectangular subset of operator entries. opSubsRef(OP,ROWIDX,COLIDX) returns subset of entries indicated by the ROWIDX and COLIDX. The indices are either integer or logical vectors, or the ':' string denoting the entire range. opSum
opSum Addition of two operators. opSum(A,B) creates a compound operator representing (A + B). opToepGauss
opToepGauss Toeplitz matrix with Gaussian entries. OP = opToepGauss(M,N,TYPE,NORMALIZED) creates an M by N Toeplitz matrix with Gaussian entries. TYPE can either be 'toeplitz' or 'circular'. For the 'toeplitz' type matrix, m+n-1 different generating entries are generated at random, whereas for the 'circular' matrix, only max(m,n) are needed. When the TYPE field is empty [], or not specified 'toeplitz' is chosen by default. Setting the NORMALIZED flag scales the columns of the Toeplitz matrix to unit norm. Multiplication is implemented using the fast Fourier transform. The use of such matrices in compressed sensing was suggested by: [1] W. U. Bajwa, J. D. Haupt, G. M. Raz, S. J. Wright, and R. D. Nowak, Toeplitz-Structured Compressed Sensing Matrices, IEEE/SP 14th Workshop on Statistical Signal Processing, pp. 294-298, August 2007. opToepSign
opToepSign Toeplitz matrix with random sign entries. OP = opToepSign(M,N,TYPE,NORMALIZED) creates an M by N Toeplitz matrix with random +1 and -1 entries. TYPE can either be 'toeplitz' or 'circular'. For the 'toeplitz' type matrix, m+n-1 different generating entries are generated at random, whereas for the 'circular' matrix, only max(m,n) are needed. When the TYPE field is empty [], or not specified 'toeplitz' is chosen by default. Setting the NORMALIZED flag scales the columns of the Toeplitz matrix to unit norm. Multiplication is implemented using the fast Fourier transform. The use of such matrices in compressed sensing was suggested by: [1] W. U. Bajwa, J. D. Haupt, G. M. Raz, S. J. Wright, and R. D. Nowak, Toeplitz-Structured Compressed Sensing Matrices, IEEE/SP 14th Workshop on Statistical Signal Processing, pp. 294-298, August 2007. opToeplitz
opToeplitz Toeplitz matrix. OP = opToeplitz(R) creates an N-by-N circular Toeplitz operator from the N-vector R. The entries of R prescribe the first row of the operator. OP = opToeplitz(C,R) creates an M-by-N Toeplitz operator where M = length(C) and N = length(R). The entries of C prescribe the first column of the operator, and likewise, R prescribes the first row. The above calls are nearly idential to Matlab's built-in TOEPLITZ function. Additionally, each call above accepts an optional logical flag that indicates if the column are scaled to have unit 2-norm length: OP = opToeplitz(R,NORMALIZED) OP = opToeplitz(C,R,NORMALIZED) Multiplication in either mode is implemented using the fast Fourier transform opTranspose
opTranspose Transpose of an operator. opTranspose(OP) returns the tranpose of OP. opUnaryMinus
opUnaryMinus Negation of an operator. opUnaryMinus(OP) returns -OP. opWavelet
opWavelet Wavelet operator. opWavelet(N) creates a Wavelet transform for 1-dimensional signals of size N. The wavelet transformation is computed using the Rice Wavelet Toolbox. opWavelet(N,FAMILY) additionally specifies the FAMILY for the wavelet. Supported values for FAMILY are 'Daubechies' and 'Haar'. opWavelet(N,FAMILY,FILTER,LEVELS,REDUNDANT,TYPE) allows for four additional parameters: FILTER (default 8) specifies the filter length, which must be even. LEVELS (default 5) gives the number of levels in the transformation. P does not need to be divisible by 2^LEVELS. However, if LEVELS is bigger than LOG2(P), then LEVELS is adjusted to be equal to FLOOR(LOG2(P)). The Boolean field REDUNDANT (default false) indicates whether the wavelet is redundant. TYPE (default 'min') indictates what type of solution is desired; 'min' for minimum phase, 'max' for maximum phase, and 'mid' for mid-phase solutions. The opWavelet operator is linear but not orthogonal. Therefore, the transpose of the operator is not the inverse operator. However, the inverse of the operator can be obtained through a left-inverse operation. For example: W = opWavelet(...) y = W*x if z = W'*y, then z ~= x but, u = W\y, then u = x opWavelet2
OPWAVELET Wavelet operator. opWavelet(P,Q,FAMILY) creates a Wavelet operator of given FAMILY for signals of size P-by-Q. The wavelet transformation is computed using the Rice Wavelet Toolbox. The values supported for FAMILY are 'Daubechies' and 'Haar'. opWavelet(P,Q,FAMILY,FILTER,LEVELS,REDUNDANT,TYPE) allows for four additional parameters: FILTER (default 8) specifies the filter length, which must be even. LEVELS (default 5) gives the number of levels in the transformation. P and Q do not need to be divisible by 2^LEVELS. However, if LEVELS is bigger than LOG2(MIN(P,Q)), then LEVELS is adjusted to be equal to FLOOR(LOG2(MIN(P,Q))). The Boolean field REDUNDANT (default false) indicates whether the wavelet is redundant. TYPE (default 'min') indictates what type of solution is desired; 'min' for minimum phase, 'max' for maximum phase, and 'mid' for mid-phase solutions. The opWavelet operator is linear but not orthogonal. Therefore, the transpose of the operator is not the inverse operator. However, the inverse of the operator can be obtained through a left-inverse operation. For example: W = opWavelet(...) y = W*x if z = W'*y, then z ~= x but, u = W\y, then u = x opWindow
opWindow Diagonal window matrix. opWindow(D) creates an operator for multiplication by the diagonal matrix with D on its diagonal. opWindow(N,FAMILY,VARARGIN) creates an N by N diagonal matrix with diagonal entries given by the window function of given family and optional additional parameters. Family Optional parameters --------------------- ---------------------------------------- Bartlett - Bartlett-Hann - Blackman alpha = 0.16; Blackman-Harris - Blackman-Nuttall - Bohman - Cauchy alpha = 3; Cos (see Cosine) Cosine alpha = 1; Dirichlet - Flattop - Gauss (see Gaussian) Gaussian alpha = 2.5; Hamming - Hann - Hann-Poisson alpha = 1; Kaiser alpha = 0.5; Kaiser-Bessel Derived alpha = 0.5; KBD (see Kaiser-Bessel Derived) Lanczos alpha = 1; Triangle - Parzen - Poisson alpha = 1; Rectangle (see Dirichlet) Sinc (see Lanczos) Tukey alpha = 0.5; Uniform (see Dirichlet) Valle-Poussin (see Parzen) Weierstrasss (see Gaussian) opZeros
opZeros Operator equivalent to zeros function. opZeros(M,N) creates an operator corresponding to an M-by-N matrix of zeros. If parameter N is omitted it is set to M. opiChol
opiChol Operator representing the incomplete Cholesky factorization of a symmetric and definitematrix with optional iterative refinement. opiChol(A) creates an operator for multiplication by an approximate inverse of the matrix A implicitly represented by its incomplete Cholesky factorization. Optionally, iterative refinement is performed. Note that A is an explicit matrix. Additional input arguments are passed directly to ichol() with the exception of the 'shape' attribute. The following attributes may be changed by the user: * nitref : the maximum number of iterative refinement steps (3) * itref_tol : iterative refinement tolerance (1.0e-8) * force_itref : force iterative refinement (false) Dominique Orban <dominique.orban@gerad.ca>, 2014. Copyright 2009, Ewout van den Berg and Michael P. Friedlander See the file COPYING.txt for full copyright information. Use the command 'spot.gpl' to locate this file. opiLU
opiLU Operator representing the incomplete LU factorization of a matrix with optional iterative refinement. If the matrix is known to be symmetric and definite, opiChol will be more efficient. opiLU(A) creates an operator for multiplication by an approximate inverse of the matrix A implicitly represented by its incomplete LU factorization. Optionally, iterative refinement is performed. Note that A is an explicit matrix. Additional input arguments are passed directly to ilu(). The following attributes may be changed by the user: * nitref : the maximum number of iterative refinement steps (3) * itref_tol : iterative refinement tolerance (1.0e-8) * force_itref : force iterative refinement (false) Dominique Orban <dominique.orban@gerad.ca>, 2014. Copyright 2009, Ewout van den Berg and Michael P. Friedlander See the file COPYING.txt for full copyright information. Use the command 'spot.gpl' to locate this file. |