User Documentation
 All Files Functions Groups
sparse_linear_systems.sql_in
1 /* ----------------------------------------------------------------------- *//**
2  *
3  * @file sparse_linear_sytems.sql_in
4  *
5  * @brief SQL functions for linear systems
6  * @date January 2011
7  *
8  * @sa Computes the solution of a consistent linear system
9  *
10  *//* ----------------------------------------------------------------------- */
11 
12 m4_include(`SQLCommon.m4')
13 
14 
15 /**
16 @addtogroup grp_sparse_linear_solver
17 
18 <div class ="toc"><b>Contents</b>
19 <ul>
20 <li class="level1"><a href="#sls_about">About</a></li>
21 <li class="level1"><a href="#sls_online_help">Online Help</a></li>
22 <li class="level1"><a href="#sls_function">Function Syntax</a></li>
23 <li class="level1"><a href="#sls_args">Arguments</a></li>
24 <li class="level1"><a href="#sls_opt_params">Optimizer Parameters</a></li>
25 <li class="level1"><a href="#sls_output">Output Tables</a></li>
26 <li class="level1"><a href="#sls_examples">Examples</a></li>
27 </ul>
28 </div>
29 
30 @anchor sls_about
31 @about
32 
33 The sparse linear systems module implements solution methods for systems of a consistent
34 linear equations. Systems of linear equations take the form:
35 \f[
36  Ax = b
37 \f]
38 
39 where \f$x \in \mathbb{R}^{n}\f$, \f$A \in \mathbb{R}^{m \times n} \f$ and \f$b \in \mathbb{R}^{m}\f$.
40 This module accepts sparse matrix input formats for \f$A\f$ and \f$b\f$.
41 We assume that there are no rows of \f$A\f$ where all elements are zero.
42 
43 \note Algorithms with fail if there is an row of the input matrix containing all zeros.
44 
45 The algorithms implemented in this module can handle large sparse
46 square linear systems. Currently, the algorithms implemented in this module
47 solve the linear system using direct or iterative methods.
48 
49 @anchor sls_online_help
50 @par Online Help
51 
52 View short help messages using the following statements:
53 @verbatim
54 -- Summary of sparse linear systems
55 SELECT madlib.linear_solver_sparse();
56 
57 -- Function syntax and output table format
58 SELECT madlib.linear_solver_sparse('usage');
59 
60 -- Syntax for direct methods
61 SELECT madlib.linear_solver_sparse('direct');
62 
63 -- Syntax for iterative methods
64 SELECT madlib.linear_solver_sparse('iterative');
65 @endverbatim
66 
67 @anchor sls_function
68 @par Syntax
69 
70 <pre>
71 linear_solver_sparse(tbl_source_lhs, tbl_source_rhs, tbl_result,
72  row_id, LHS, RHS, grouping_cols := NULL,
73  optimizer := 'direct',
74  optimizer_params := 'algorithm = llt');
75 </pre>
76 
77 @anchor sls_args
78 @par Arguments
79 
80 
81 <DL class="arglist">
82 <DT>LHS</DT>
83 <DD>Text value. The name of the table containing the left hand side matrix.
84 For the LHS matrix, the input data is expected to be of the following form:
85 <pre>{TABLE|VIEW} <em>sourceName</em> (
86  ...
87  <em>row_id</em> FLOAT8,
88  <em>col_id</em> FLOAT8,
89  <em>value</em> FLOAT8,
90  ...
91 )</pre>
92 
93 
94 Here, each row represents a single equation using. The <em> rhs </em> refers
95 to the right hand side of the equations while the <em> lhs</em>
96 refers to the multipliers on the variables on the left hand side of the same
97 equations.
98 </DD>
99 
100 <DT>RHS</DT>
101 <DD>Text value. The name of the table containing the right hand side vector.
102 For the RHS matrix, the input data is expected to be of the following form:
103 <pre>{TABLE|VIEW} <em>sourceName</em> (
104  ...
105  <em>row_id</em> FLOAT8,
106  <em>value</em> FLOAT8
107  ...
108 )</pre>
109 Here, each row represents a single equation using. The <em> rhs </em> refers
110 to the right hand side of the equations while the <em> lhs</em>
111 refers to the multipliers on the variables on the left hand side of the same
112 equations.
113 
114 Output is stored in the <em>tbl_result</em>
115 @verbatim
116  tbl_result | Data Types
117 --------------------|---------------------
118 solution | DOUBLE PRECISION[]
119 residual_norm | DOUBLE PRECISIOn
120 iters | INTEGER
121 @endverbatim
122 
123 @par Syntax
124 
125 <pre>
126 SELECT sparse_linear_sytems('tbl_source_lhs', 'tbl_source_rhs', 'tbl_result',
127  'row_id', 'LHS', 'RHS', NULL, 'direct', 'algorithm = householderqr');
128 </pre>
129 
130 <DL>
131 <DT>tbl_source_lhs</DT>
132 <DD>Text value. The name of the table containing the LHS matrix.</DD>
133 
134 <DT>tbl_source_rhs</DT>
135 <DD>Text value. The name of the table containing the RHS vector.</DD>
136 
137 <DT>tbl_result</DT>
138 <DD>Text value. The name of the table where the output is saved.</DD>
139 
140 <DT>lhs_row_id</DT>
141 <DD>Text value. The name of the column storing the 'row id' of the equations.</DD>
142 
143 \note For a system with N equations, the row_id's must be a continuous
144 range of integers from \f$ 0 \ldots n-1 \f$.
145 
146 
147 <DT>lhs_col_id</DT>
148 <DD>Text value. The name of the column (in tbl_source_lhs) storing the 'col id' of the equations.</DD>
149 
150 <DT>lhs_value</DT>
151 <DD>Text value. The name of the column (in tbl_source_lhs) storing the 'value' of the equations.</DD>
152 
153 
154 <DT>rhs_row_id</DT>
155 <DD>Text value. The name of the column (in tbl_source_rhs) storing the 'col id' of the equations.</DD>
156 
157 <DT>rhs_value</DT>
158 <DD>Text value. The name of the column (in tbl_source_rhs) storing the 'value' of the equations.</DD>
159 
160 <DT>num_vars</DT>
161 <DD>Integer value. The number of variables in the linear system.
162 equations.</DD>
163 
164 <DT>grouping_col (optional) </DT>
165 <DD>Text value. Group by columns. Default: NULL.</DD>
166 \note The grouping columns is a placeholder in MADlib V1.1
167 
168 <DT>optimizer (optional) </DT>
169 <DD>Text value. Type of optimizer. Default: 'direct'.</DD>
170 
171 <DT>optimizer_params (optional) </DT>
172 <DD>Text value. Optimizer specific parameters. Default: NULL.</DD>
173 </DL>
174 
175 @anchor sls_opt_params
176 @par Optimizer Parameters
177 
178 For each optimizer, there are specific parameters that can be tuned
179 for better performance.
180 
181 <DL class="arglist">
182 <DT>algorithm (default: ldlt)</dT>
183 <DD>
184 
185  There are several algorithms that can be classified as 'direct' methods
186  of solving linear systems. Madlib functions provide various algorithmic
187  options available for users.
188 
189  The following table provides a guideline on the choice of algorithm based
190  on conditions on the A matrix, speed of the algorithms and numerical stability.
191 
192  Algorithm | Contitions on A | Speed | Memory
193  ----------------------------------------------------------
194  llt | Sym. Pos Def | ++ | ---
195  ldlt | Sym. Pos Def | ++ | ---
196 
197  For speed '++' is faster than '+', which is faster than '-'.
198  For accuracy '+++' is better than '++'.
199  For memory, '-' uses less memory than '--'.
200 
201  Note: ldlt is often preferred over llt
202 
203 
204  There are several algorithms that can be classified as 'iterative' methods
205  of solving linear systems. Madlib functions provide various algorithmic
206  options available for users.
207 
208  The following table provides a guideline on the choice of algorithm based
209  on conditions on the A matrix, speed of the algorithms and numerical stability.
210 
211 
212  Algorithm | Contitions on A | Speed | Memory | Convergence
213  ----------------------------------------------------------------------
214  cg-mem | Sym. Pos Def | +++ | - | ++
215  bicgstab-mem | Square | ++ | - | +
216  precond-cg-mem | Sym. Pos Def | ++ | - | +++
217  precond-bicgstab-mem | Square | + | - | ++
218 
219  For memory, '-' uses less memory than '--'.
220  For speed, '++' is faster than '+'.
221 
222 Details:
223 -# <b>cg-mem: </b> In memory conjugate gradient with diagonal preconditioners.
224 -# <b>bicgstab-mem: </b> Bi-conjugate gradient (equivalent to performing CG on the least squares formulation of Ax=b) with incomplete LU preconditioners.
225 -# <b>precond-cg-mem:</b> In memory conjugate gradient with diagonal preconditioners.
226 -# <b>bicgstab-mem: </b> Bi-conjugate gradient (equivalent to performing CG on the least squares formulation of Ax=b) with incomplete LU preconditioners.
227 </DD>
228 
229 <DT> toler (default: 1e-5)</dT>
230 <DD> Termination tolerance (applicable only for iterative methods) which
231 determines the stopping criterion (with respect to residual norm) for iterative methods.
232 </DD>
233 
234 </DL>
235 
236 
237 @anchor sls_output
238 @par Output statistics
239 
240 <DL class="arglist">
241 <DT>solution</dT>
242 <DD>
243 The solution is an array (of double precision) with the variables in the same
244 order as that provided as input in the 'left_hand_side' column name of the
245 'source_table'
246 </DD>
247 
248 <DT>residual_norm</dT>
249 <DD>
250 Computes the scaled residual norm, defined as \f$ \frac{|Ax - b|}{|b|} \f$
251 gives the user an indication of the accuracy of the solution.
252 </DD>
253 
254 <DT>iters</dT>
255 <DD>`
256 The number of iterations required by the algorithm (only applicable for
257  iterative algorithms) . The output is NULL for 'direct' methods.
258 </DD>
259 
260 </DL>
261 
262 
263 @anchor sls_examples
264 @examp
265 
266 -# Create the sample data set:
267 \verbatim
268 sql> DROP TABLE IF EXISTS sparse_linear_systems_lhs;
269  CREATE TABLE sparse_linear_systems_lhs (
270  rid INTEGER NOT NULL,
271  cid INTEGER,
272  val DOUBLE PRECISION
273  );
274 
275  DROP TABLE IF EXISTS sparse_linear_systems_rhs;
276  CREATE TABLE sparse_linear_systems_rhs (
277  rid INTEGER NOT NULL,
278  val DOUBLE PRECISION
279  );
280 
281 
282  INSERT INTO sparse_linear_systems_lhs(rid, cid, val) VALUES
283  (0, 0, 1),
284  (1, 1, 1),
285  (2, 2, 1),
286  (3, 3, 1);
287 
288  INSERT INTO sparse_linear_systems_rhs(rid, val) VALUES
289  (0, 10),
290  (1, 20),
291  (2, 30);
292 \endverbatim
293 
294 -# Solve the linear systems with default parameters
295 \verbatim
296 sql> SELECT madlib.linear_solver_sparse(
297  'sparse_linear_systems_lhs',
298  'sparse_linear_systems_rhs',
299  'output_table',
300  'rid',
301  'cid',
302  'val',
303  'rid',
304  'val',
305  4);
306 \endverbatim
307 
308 -# Obtain the output from the output table
309 \verbatim
310 sql> SELECT * FROM output_table;
311 --------------------+-------------------------------------
312 solution | {10,20,30,0}
313 residual_norm | 0
314 iters | NULL
315 \endverbatim
316 
317 
318 -# Chose a different algorithm than the default algorithm
319 \verbatim
320 drop table if exists result_table;
321 
322 SELECT madlib.linear_solver_sparse(
323  'sparse_linear_systems_lhs',
324  'sparse_linear_systems_rhs',
325  'output_table',
326  'rid',
327  'cid',
328  'val',
329  'rid',
330  'val',
331  4,
332  NULL,
333  'direct',
334  'algorithm=llt');
335 
336 -# Chose a different algorithm than the default algorithm
337 \verbatim
338 drop table if exists result_table;
339 madlib.linear_solver_sparse(
340  'sparse_linear_systems_lhs',
341  'sparse_linear_systems_rhs',
342  'output_table',
343  'rid',
344  'cid',
345  'val',
346  'rid',
347  'val',
348  4,
349  NULL,
350  'iterative',
351  'algorithm=cg-mem, toler=1e-5');
352 \endverbatim
353 
354 @sa File sparse_linear_sytems.sql_in documenting the SQL functions.
355 
356 @internal
357 @sa Namespace \ref madlib::modules::linear_systems documenting the implementation in C++
358 @endinternal
359 */
360 
361 
362 ------------------ Linear Systems ------------------------------
363 
364 CREATE TYPE MADLIB_SCHEMA.sparse_linear_solver_result AS (
365  solution DOUBLE PRECISION[],
366  residual_norm DOUBLE PRECISION,
367  iters INTEGER
368 );
369 
370 ------------------ In memory Iterative ------------------------------
371 
372 CREATE OR REPLACE FUNCTION MADLIB_SCHEMA.sparse_inmem_iterative_linear_system_transition(
373  state DOUBLE PRECISION[],
374  row_id INTEGER,
375  col_id INTEGER,
376  value DOUBLE PRECISION,
377  b DOUBLE PRECISION,
378  num_eqs INTEGER,
379  num_vars INTEGER,
380  nnz INTEGER,
381  algorithm INTEGER,
382  maxIter INTEGER,
383  termToler DOUBLE PRECISION)
384 RETURNS DOUBLE PRECISION[]
385 AS 'MODULE_PATHNAME'
386 LANGUAGE C IMMUTABLE STRICT;
387 
388 CREATE OR REPLACE FUNCTION MADLIB_SCHEMA.sparse_inmem_iterative_linear_system_merge_states(
389  state1 DOUBLE PRECISION[],
390  state2 DOUBLE PRECISION[])
391 RETURNS DOUBLE PRECISION[]
392 AS 'MODULE_PATHNAME'
393 LANGUAGE C IMMUTABLE STRICT;
394 
395 
396 CREATE OR REPLACE FUNCTION MADLIB_SCHEMA.sparse_inmem_iterative_linear_system_final(
397  state DOUBLE PRECISION[])
398 RETURNS MADLIB_SCHEMA.sparse_linear_solver_result
399 AS 'MODULE_PATHNAME'
400 LANGUAGE C IMMUTABLE STRICT;
401 
402 
403 /**
404  * @brief Solve a system of linear equations using the inmem_iterative method
405  *
406  * @param row_id Column containing the row_id
407  * @param col_id Column containing the col_id
408  * @param value Value of the LHS matrix
409  * @param right_hand_side Column containing the right hand side of the system
410  * @param numEquations Number of equations
411  * @param algorithm Algorithm used for the sparse linear solver
412  * @param maxIter Maximum number of iterations
413  * @param termToler Termination tolerance
414  *
415  *
416  * @return A composite value:
417  * - <tt>solution FLOAT8[] </tt> - Array of marginal effects
418  * - <tt>residual_norm FLOAT8</tt> - Norm of the residual
419  * - <tt>iters INTEGER</tt> - Iterations taken
420  *
421  * @usage
422  * - Get all the diagnostic statistics:\n
423  *
424  * <pre> SELECT linear_system_sparse(<em>row_id</em>,
425  * <em>left_hand_side</em>,
426  * <em> right_hand_side </em>,
427  * <em> numEquations </em>)
428  * FROM <em>dataTable</em>;
429  * </pre>
430  */
431 
432 CREATE AGGREGATE MADLIB_SCHEMA.sparse_inmem_iterative_linear_system(
433  /*+ "row_id" */ INTEGER,
434  /*+ "col_id" */ INTEGER,
435  /*+ "value" */ DOUBLE PRECISION,
436  /*+ "right_hand_side" */ DOUBLE PRECISION,
437  /*+ "numEquations" */ INTEGER,
438  /*+ "numVars" */ INTEGER,
439  /*+ "nnz" */ INTEGER,
440  /*+ "algorithm" */ INTEGER,
441  /*+ "maxIter" */ INTEGER,
442  /*+ "termToler" */ DOUBLE PRECISION)(
443  STYPE=DOUBLE PRECISION[],
444  SFUNC=MADLIB_SCHEMA.sparse_inmem_iterative_linear_system_transition,
445  m4_ifdef(`__GREENPLUM__',`PREFUNC=MADLIB_SCHEMA.sparse_inmem_iterative_linear_system_merge_states,')
446  FINALFUNC=MADLIB_SCHEMA.sparse_inmem_iterative_linear_system_final,
447  INITCOND='{0,0,0,0,0,0,0,0}'
448 );
449 
450 ------------------ Direct Method ------------------------------
451 
452 CREATE OR REPLACE FUNCTION MADLIB_SCHEMA.sparse_direct_linear_system_transition(
453  state DOUBLE PRECISION[],
454  row_id INTEGER,
455  col_id INTEGER,
456  value DOUBLE PRECISION,
457  b DOUBLE PRECISION,
458  num_eqs INTEGER,
459  num_vars INTEGER,
460  nnz INTEGER,
461  algorithm INTEGER)
462 RETURNS DOUBLE PRECISION[]
463 AS 'MODULE_PATHNAME'
464 LANGUAGE C IMMUTABLE STRICT;
465 
466 CREATE OR REPLACE FUNCTION MADLIB_SCHEMA.sparse_direct_linear_system_merge_states(
467  state1 DOUBLE PRECISION[],
468  state2 DOUBLE PRECISION[])
469 RETURNS DOUBLE PRECISION[]
470 AS 'MODULE_PATHNAME'
471 LANGUAGE C IMMUTABLE STRICT;
472 
473 
474 CREATE OR REPLACE FUNCTION MADLIB_SCHEMA.sparse_direct_linear_system_final(
475  state DOUBLE PRECISION[])
476 RETURNS MADLIB_SCHEMA.sparse_linear_solver_result
477 AS 'MODULE_PATHNAME'
478 LANGUAGE C IMMUTABLE STRICT;
479 
480 
481 /**
482  * @brief Solve a system of linear equations using the direct method
483  *
484  * @param row_id Column containing the row_id
485  * @param col_id Column containing the col_id
486  * @param value Value of the LHS matrix
487  * @param right_hand_side Column containing the right hand side of the system
488  * @param numEquations Number of equations
489  * @param algorithm Algorithm used for the sparse linear solver
490  *
491  *
492  * @return A composite value:
493  * - <tt>solution FLOAT8[] </tt> - Array of marginal effects
494  * - <tt>residual_norm FLOAT8</tt> - Norm of the residual
495  * - <tt>iters INTEGER</tt> - Iterations taken
496  *
497  * @usage
498  * - Get all the diagnostic statistics:\n
499  *
500  * <pre> SELECT linear_system_sparse(<em>row_id</em>,
501  * <em>left_hand_side</em>,
502  * <em> right_hand_side </em>,
503  * <em> numEquations </em>)
504  * FROM <em>dataTable</em>;
505  * </pre>
506  */
507 
508 CREATE AGGREGATE MADLIB_SCHEMA.sparse_direct_linear_system(
509  /*+ "row_id" */ INTEGER,
510  /*+ "col_id" */ INTEGER,
511  /*+ "value" */ DOUBLE PRECISION,
512  /*+ "right_hand_side" */ DOUBLE PRECISION,
513  /*+ "numEquations" */ INTEGER,
514  /*+ "numVars" */ INTEGER,
515  /*+ "nnz" */ INTEGER,
516  /*+ "algorithm" */ INTEGER)(
517  STYPE=DOUBLE PRECISION[],
518  SFUNC=MADLIB_SCHEMA.sparse_direct_linear_system_transition,
519  m4_ifdef(`__GREENPLUM__',`PREFUNC=MADLIB_SCHEMA.sparse_direct_linear_system_merge_states,')
520  FINALFUNC=MADLIB_SCHEMA.sparse_direct_linear_system_final,
521  INITCOND='{0,0,0,0,0,0}'
522 );
523 
524 
525 --------------------------- Interface ----------------------------------
526 
527 /**
528  * @brief Help function, to print out the supported families
529  */
530 CREATE OR REPLACE FUNCTION MADLIB_SCHEMA.linear_solver_sparse(
531  input_string VARCHAR
532 )
533 RETURNS VARCHAR AS $$
534 PythonFunction(linear_systems, sparse_linear_systems, linear_solver_sparse_help)
535 $$ LANGUAGE plpythonu;
536 
537 CREATE OR REPLACE FUNCTION MADLIB_SCHEMA.linear_solver_sparse()
538 RETURNS VARCHAR AS $$
539 BEGIN
540  RETURN MADLIB_SCHEMA.linear_solver_sparse('');
541 END;
542 $$ LANGUAGE plpgsql VOLATILE;
543 
544 
545 /**
546  @brief A wrapper function for the various marginal linear_systemsion analyzes.
547  *
548  * @param lhs_table String identifying the input A matrix
549  * @param rhs_table String identifying the input b vector
550  * @param out_table String identifying the output table to be created
551  * @param lhs_row_id Column name containing the LHS of the equations
552  * @param lhs_col_id Column name containing the LHS of the equations
553  * @param lhs_value Column name containing the LHS of the equations
554  * @param rhs_row_id Column name containing the RHS of the equations
555  * @param rhs_value Column name containing the RHS of the equations
556  * @param num_vars Number of variables in the system
557  * @param grouping_sols Columns to group the linear systems by
558  * @param optimzer Optimizer to be used
559  * @param optimzer_options Optimizer options used
560  *
561  *
562  * @return void
563  *
564  * @usage
565  * For function summary information. Run
566  * sql> select linear_solver_sparse('help');
567  * OR
568  * sql> select linear_solver_sparse();
569  * OR
570  * sql> select linear_solver_sparse('?');
571  * For function usage information. Run
572  * sql> select linear_solver_sparse('usage');
573  *
574  */
575 
576 CREATE OR REPLACE FUNCTION MADLIB_SCHEMA.linear_solver_sparse(
577  lhs_table VARCHAR -- name of input lhs table
578  , rhs_table VARCHAR -- name of input rhs table
579  , out_table VARCHAR -- name of output table
580  , lhs_row_id VARCHAR -- column name with row_id
581  , lhs_col_id VARCHAR -- column name with col_id
582  , lhs_value VARCHAR -- column name with value
583  , rhs_row_id VARCHAR -- rhs row_id
584  , rhs_value VARCHAR -- rhs value
585  , num_vars INTEGER -- Number of variables
586  , grouping_cols VARCHAR -- name of columns to group by
587  , optimizer VARCHAR -- Name of the optimizer
588  , optimizer_options VARCHAR -- Optimal parameters of the optimizer
589  )
590 RETURNS VOID AS $$
591 PythonFunction(linear_systems, sparse_linear_systems, linear_solver_sparse)
592 $$ LANGUAGE plpythonu;
593 
594 
595 
596 -- Default Variable calls for linear_solver_sparse
597 ------------------------------------------------------------------------------
598 
599 /**
600  * @brief Marginal effects with default variables
601  **/
602 CREATE OR REPLACE FUNCTION MADLIB_SCHEMA.linear_solver_sparse(
603  lhs_table VARCHAR -- name of input lhs table
604  , rhs_table VARCHAR -- name of input rhs table
605  , out_table VARCHAR -- name of output table
606  , lhs_row_id VARCHAR -- column name with row_id
607  , lhs_col_id VARCHAR -- column name with col_id
608  , lhs_value VARCHAR -- column name with value
609  , rhs_row_id VARCHAR -- rhs row_id
610  , rhs_value VARCHAR -- rhs value
611  , num_vars INTEGER -- Number of variables
612  )
613 RETURNS VOID AS $$
614 BEGIN
615  PERFORM MADLIB_SCHEMA.linear_solver_sparse(
616  lhs_table,
617  rhs_table,
618  out_table,
619  lhs_row_id,
620  lhs_col_id,
621  lhs_value,
622  rhs_row_id,
623  rhs_value,
624  num_vars,
625  NULL,
626  'direct',
627  NULL);
628 
629 END;
630 $$ LANGUAGE plpgsql VOLATILE;
631 
632 
633 /**
634  * @brief Marginal effects with default variables
635  **/
636 CREATE OR REPLACE FUNCTION MADLIB_SCHEMA.linear_solver_sparse(
637  lhs_table VARCHAR -- name of input lhs table
638  , rhs_table VARCHAR -- name of input rhs table
639  , out_table VARCHAR -- name of output table
640  , lhs_row_id VARCHAR -- column name with row_id
641  , lhs_col_id VARCHAR -- column name with col_id
642  , lhs_value VARCHAR -- column name with value
643  , rhs_row_id VARCHAR -- rhs row_id
644  , rhs_value VARCHAR -- rhs value
645  , num_vars INTEGER -- Number of variables
646  , grouping_cols VARCHAR -- name of columns to group by
647 )
648 RETURNS VOID AS $$
649 BEGIN
650  PERFORM MADLIB_SCHEMA.linear_solver_sparse(
651  lhs_table,
652  rhs_table,
653  out_table,
654  lhs_row_id,
655  lhs_col_id,
656  lhs_value,
657  rhs_row_id,
658  rhs_value,
659  num_vars,
660  grouping_cols,
661  'direct',
662  NULL);
663 END;
664 $$ LANGUAGE plpgsql VOLATILE;
665 
666 /**
667  * @brief Marginal effects with default variables
668  **/
669 CREATE OR REPLACE FUNCTION MADLIB_SCHEMA.linear_solver_sparse(
670  lhs_table VARCHAR -- name of input lhs table
671  , rhs_table VARCHAR -- name of input rhs table
672  , out_table VARCHAR -- name of output table
673  , lhs_row_id VARCHAR -- column name with row_id
674  , lhs_col_id VARCHAR -- column name with col_id
675  , lhs_value VARCHAR -- column name with value
676  , rhs_row_id VARCHAR -- rhs row_id
677  , rhs_value VARCHAR -- rhs value
678  , num_vars INTEGER -- Number of variables
679  , grouping_cols VARCHAR -- name of columns to group by
680  , optimizer VARCHAR -- Name of the optimizer
681  )
682 RETURNS VOID AS $$
683 BEGIN
684  PERFORM MADLIB_SCHEMA.linear_solver_sparse(
685  lhs_table,
686  rhs_table,
687  out_table,
688  lhs_row_id,
689  lhs_col_id,
690  lhs_value,
691  rhs_row_id,
692  rhs_value,
693  num_vars,
694  grouping_cols,
695  optimizer,
696  NULL);
697 END;
698 $$ LANGUAGE plpgsql VOLATILE;