1use crate::enums;
165use crate::Value;
166use ffi::FFI;
167
168use types::complex::FFFI;
169
170#[doc(alias = "gsl_linalg_LU_decomp")]
186pub fn LU_decomp(
187 a: &mut crate::MatrixF64,
188 p: &mut ::Permutation,
189 signum: &mut i32,
190) -> Result<(), Value> {
191 let ret = unsafe { sys::gsl_linalg_LU_decomp(a.unwrap_unique(), p.unwrap_unique(), signum) };
192 result_handler!(ret, ())
193}
194
195#[doc(alias = "gsl_linalg_complex_LU_decomp")]
211pub fn complex_LU_decomp(
212 a: &mut crate::MatrixComplexF64,
213 p: &mut ::Permutation,
214 signum: &mut i32,
215) -> Result<(), Value> {
216 let ret =
217 unsafe { sys::gsl_linalg_complex_LU_decomp(a.unwrap_unique(), p.unwrap_unique(), signum) };
218 result_handler!(ret, ())
219}
220
221#[doc(alias = "gsl_linalg_LU_solve")]
223pub fn LU_solve(
224 lu: &::MatrixF64,
225 p: &::Permutation,
226 b: &::VectorF64,
227 x: &mut crate::VectorF64,
228) -> Result<(), Value> {
229 let ret = unsafe {
230 sys::gsl_linalg_LU_solve(
231 lu.unwrap_shared(),
232 p.unwrap_shared(),
233 b.unwrap_shared(),
234 x.unwrap_unique(),
235 )
236 };
237 result_handler!(ret, ())
238}
239
240#[doc(alias = "gsl_linalg_complex_LU_solve")]
242pub fn complex_LU_solve(
243 lu: &::MatrixComplexF64,
244 p: &::Permutation,
245 b: &::VectorComplexF64,
246 x: &mut crate::VectorComplexF64,
247) -> Result<(), Value> {
248 let ret = unsafe {
249 sys::gsl_linalg_complex_LU_solve(
250 lu.unwrap_shared(),
251 p.unwrap_shared(),
252 b.unwrap_shared(),
253 x.unwrap_unique(),
254 )
255 };
256 result_handler!(ret, ())
257}
258
259#[doc(alias = "gsl_linalg_LU_svx")]
262pub fn LU_svx(lu: &::MatrixF64, p: &::Permutation, x: &mut crate::VectorF64) -> Result<(), Value> {
263 let ret =
264 unsafe { sys::gsl_linalg_LU_svx(lu.unwrap_shared(), p.unwrap_shared(), x.unwrap_unique()) };
265 result_handler!(ret, ())
266}
267
268#[doc(alias = "gsl_linalg_complex_LU_svx")]
271pub fn complex_LU_svx(
272 lu: &::MatrixComplexF64,
273 p: &::Permutation,
274 x: &mut crate::VectorComplexF64,
275) -> Result<(), Value> {
276 let ret = unsafe {
277 sys::gsl_linalg_complex_LU_svx(lu.unwrap_shared(), p.unwrap_shared(), x.unwrap_unique())
278 };
279 result_handler!(ret, ())
280}
281
282#[doc(alias = "gsl_linalg_LU_refine")]
285pub fn LU_refine(
286 a: &::MatrixF64,
287 lu: &::MatrixF64,
288 p: &::Permutation,
289 b: &::VectorF64,
290 x: &mut crate::VectorF64,
291 residual: &mut crate::VectorF64,
292) -> Result<(), Value> {
293 let ret = unsafe {
294 sys::gsl_linalg_LU_refine(
295 a.unwrap_shared(),
296 lu.unwrap_shared(),
297 p.unwrap_shared(),
298 b.unwrap_shared(),
299 x.unwrap_unique(),
300 residual.unwrap_unique(),
301 )
302 };
303 result_handler!(ret, ())
304}
305
306#[doc(alias = "gsl_linalg_complex_LU_refine")]
309pub fn complex_LU_refine(
310 a: &mut crate::MatrixComplexF64,
311 lu: &::MatrixComplexF64,
312 p: &::Permutation,
313 b: &::VectorComplexF64,
314 x: &mut crate::VectorComplexF64,
315 residual: &mut crate::VectorComplexF64,
316) -> Result<(), Value> {
317 let ret = unsafe {
318 sys::gsl_linalg_complex_LU_refine(
319 a.unwrap_unique(),
320 lu.unwrap_shared(),
321 p.unwrap_shared(),
322 b.unwrap_shared(),
323 x.unwrap_unique(),
324 residual.unwrap_unique(),
325 )
326 };
327 result_handler!(ret, ())
328}
329
330#[doc(alias = "gsl_linalg_LU_invert")]
335pub fn LU_invert(
336 lu: &::MatrixF64,
337 p: &::Permutation,
338 inverse: &mut crate::MatrixF64,
339) -> Result<(), Value> {
340 let ret = unsafe {
341 sys::gsl_linalg_LU_invert(
342 lu.unwrap_shared(),
343 p.unwrap_shared(),
344 inverse.unwrap_unique(),
345 )
346 };
347 result_handler!(ret, ())
348}
349
350#[doc(alias = "gsl_linalg_complex_LU_invert")]
355pub fn complex_LU_invert(
356 lu: &::MatrixComplexF64,
357 p: &::Permutation,
358 inverse: &mut crate::MatrixComplexF64,
359) -> Result<(), Value> {
360 let ret = unsafe {
361 sys::gsl_linalg_complex_LU_invert(
362 lu.unwrap_shared(),
363 p.unwrap_shared(),
364 inverse.unwrap_unique(),
365 )
366 };
367 result_handler!(ret, ())
368}
369
370#[doc(alias = "gsl_linalg_LU_det")]
373pub fn LU_det(lu: &mut crate::MatrixF64, signum: i32) -> f64 {
374 unsafe { sys::gsl_linalg_LU_det(lu.unwrap_unique(), signum) }
375}
376
377#[doc(alias = "gsl_linalg_complex_LU_det")]
380pub fn complex_LU_det(lu: &mut crate::MatrixComplexF64, signum: i32) -> crate::ComplexF64 {
381 unsafe { sys::gsl_linalg_complex_LU_det(lu.unwrap_unique(), signum).wrap() }
382}
383
384#[doc(alias = "gsl_linalg_LU_lndet")]
387pub fn LU_lndet(lu: &mut crate::MatrixF64) -> f64 {
388 unsafe { sys::gsl_linalg_LU_lndet(lu.unwrap_unique()) }
389}
390
391#[doc(alias = "gsl_linalg_complex_LU_lndet")]
393pub fn complex_LU_lndet(lu: &mut crate::MatrixComplexF64) -> f64 {
394 unsafe { sys::gsl_linalg_complex_LU_lndet(lu.unwrap_unique()) }
395}
396
397#[doc(alias = "gsl_linalg_LU_sgndet")]
399pub fn LU_sgndet(lu: &mut crate::MatrixF64, signum: i32) -> i32 {
400 unsafe { sys::gsl_linalg_LU_sgndet(lu.unwrap_unique(), signum) }
401}
402
403#[doc(alias = "gsl_linalg_complex_LU_sgndet")]
405pub fn complex_LU_sgndet(lu: &mut crate::MatrixComplexF64, signum: i32) -> crate::ComplexF64 {
406 unsafe { sys::gsl_linalg_complex_LU_sgndet(lu.unwrap_unique(), signum).wrap() }
407}
408
409#[doc(alias = "gsl_linalg_QR_decomp")]
417pub fn QR_decomp(a: &mut crate::MatrixF64, tau: &mut crate::VectorF64) -> Result<(), Value> {
418 let ret = unsafe { sys::gsl_linalg_QR_decomp(a.unwrap_unique(), tau.unwrap_unique()) };
419 result_handler!(ret, ())
420}
421
422#[doc(alias = "gsl_linalg_QR_solve")]
425pub fn QR_solve(
426 qr: &::MatrixF64,
427 tau: &::VectorF64,
428 b: &::VectorF64,
429 x: &mut crate::VectorF64,
430) -> Result<(), Value> {
431 let ret = unsafe {
432 sys::gsl_linalg_QR_solve(
433 qr.unwrap_shared(),
434 tau.unwrap_shared(),
435 b.unwrap_shared(),
436 x.unwrap_unique(),
437 )
438 };
439 result_handler!(ret, ())
440}
441
442#[doc(alias = "gsl_linalg_QR_svx")]
445pub fn QR_svx(qr: &::MatrixF64, tau: &::VectorF64, x: &mut crate::VectorF64) -> Result<(), Value> {
446 let ret = unsafe {
447 sys::gsl_linalg_QR_svx(qr.unwrap_shared(), tau.unwrap_shared(), x.unwrap_unique())
448 };
449 result_handler!(ret, ())
450}
451
452#[doc(alias = "gsl_linalg_QR_lssolve")]
457pub fn QR_lssolve(
458 qr: &::MatrixF64,
459 tau: &::VectorF64,
460 b: &::VectorF64,
461 x: &mut crate::VectorF64,
462 residual: &mut crate::VectorF64,
463) -> Result<(), Value> {
464 let ret = unsafe {
465 sys::gsl_linalg_QR_lssolve(
466 qr.unwrap_shared(),
467 tau.unwrap_shared(),
468 b.unwrap_shared(),
469 x.unwrap_unique(),
470 residual.unwrap_unique(),
471 )
472 };
473 result_handler!(ret, ())
474}
475
476#[doc(alias = "gsl_linalg_QR_QTvec")]
479pub fn QR_QTvec(
480 qr: &::MatrixF64,
481 tau: &::VectorF64,
482 v: &mut crate::VectorF64,
483) -> Result<(), Value> {
484 let ret = unsafe {
485 sys::gsl_linalg_QR_QTvec(qr.unwrap_shared(), tau.unwrap_shared(), v.unwrap_unique())
486 };
487 result_handler!(ret, ())
488}
489
490#[doc(alias = "gsl_linalg_QR_Qvec")]
493pub fn QR_Qvec(qr: &::MatrixF64, tau: &::VectorF64, v: &mut crate::VectorF64) -> Result<(), Value> {
494 let ret = unsafe {
495 sys::gsl_linalg_QR_Qvec(qr.unwrap_shared(), tau.unwrap_shared(), v.unwrap_unique())
496 };
497 result_handler!(ret, ())
498}
499
500#[doc(alias = "gsl_linalg_QR_QTmat")]
503pub fn QR_QTmat(
504 qr: &::MatrixF64,
505 tau: &::VectorF64,
506 v: &mut crate::MatrixF64,
507) -> Result<(), Value> {
508 let ret = unsafe {
509 sys::gsl_linalg_QR_QTmat(qr.unwrap_shared(), tau.unwrap_shared(), v.unwrap_unique())
510 };
511 result_handler!(ret, ())
512}
513
514#[doc(alias = "gsl_linalg_QR_Rsolve")]
517pub fn QR_Rsolve(qr: &::MatrixF64, b: &::VectorF64, x: &mut crate::VectorF64) -> Result<(), Value> {
518 let ret = unsafe {
519 sys::gsl_linalg_QR_Rsolve(qr.unwrap_shared(), b.unwrap_shared(), x.unwrap_unique())
520 };
521 result_handler!(ret, ())
522}
523
524#[doc(alias = "gsl_linalg_QR_Rsvx")]
527pub fn QR_Rsvx(qr: &::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
528 let ret = unsafe { sys::gsl_linalg_QR_Rsvx(qr.unwrap_shared(), x.unwrap_unique()) };
529 result_handler!(ret, ())
530}
531
532#[doc(alias = "gsl_linalg_QR_unpack")]
534pub fn QR_unpack(
535 qr: &::MatrixF64,
536 tau: &::VectorF64,
537 q: &mut crate::MatrixF64,
538 r: &mut crate::MatrixF64,
539) -> Result<(), Value> {
540 let ret = unsafe {
541 sys::gsl_linalg_QR_unpack(
542 qr.unwrap_shared(),
543 tau.unwrap_shared(),
544 q.unwrap_unique(),
545 r.unwrap_unique(),
546 )
547 };
548 result_handler!(ret, ())
549}
550
551#[doc(alias = "gsl_linalg_QR_QRsolve")]
554pub fn QR_QRsolve(
555 q: &mut crate::MatrixF64,
556 r: &mut crate::MatrixF64,
557 b: &::VectorF64,
558 x: &mut crate::VectorF64,
559) -> Result<(), Value> {
560 let ret = unsafe {
561 sys::gsl_linalg_QR_QRsolve(
562 q.unwrap_unique(),
563 r.unwrap_unique(),
564 b.unwrap_shared(),
565 x.unwrap_unique(),
566 )
567 };
568 result_handler!(ret, ())
569}
570
571#[doc(alias = "gsl_linalg_QR_update")]
574pub fn QR_update(
575 q: &mut crate::MatrixF64,
576 r: &mut crate::MatrixF64,
577 mut w: crate::VectorF64,
578 v: &::VectorF64,
579) -> Result<(), Value> {
580 let ret = unsafe {
581 sys::gsl_linalg_QR_update(
582 q.unwrap_unique(),
583 r.unwrap_unique(),
584 w.unwrap_unique(),
585 v.unwrap_shared(),
586 )
587 };
588 result_handler!(ret, ())
589}
590
591#[doc(alias = "gsl_linalg_R_solve")]
593pub fn R_solve(r: &::MatrixF64, b: &::VectorF64, x: &mut crate::VectorF64) -> Result<(), Value> {
594 let ret =
595 unsafe { sys::gsl_linalg_R_solve(r.unwrap_shared(), b.unwrap_shared(), x.unwrap_unique()) };
596 result_handler!(ret, ())
597}
598
599#[doc(alias = "gsl_linalg_R_svx")]
602pub fn R_svx(r: &::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
603 let ret = unsafe { sys::gsl_linalg_R_svx(r.unwrap_shared(), x.unwrap_unique()) };
604 result_handler!(ret, ())
605}
606
607#[doc(alias = "gsl_linalg_QRPT_decomp")]
617pub fn QRPT_decomp(
618 a: &mut crate::MatrixF64,
619 tau: &mut crate::VectorF64,
620 p: &mut ::Permutation,
621 signum: &mut i32,
622 norm: &mut crate::VectorF64,
623) -> Result<(), Value> {
624 let ret = unsafe {
625 sys::gsl_linalg_QRPT_decomp(
626 a.unwrap_unique(),
627 tau.unwrap_unique(),
628 p.unwrap_unique(),
629 signum,
630 norm.unwrap_unique(),
631 )
632 };
633 result_handler!(ret, ())
634}
635
636#[doc(alias = "gsl_linalg_QRPT_decomp2")]
639pub fn QRPT_decomp2(
640 a: &::MatrixF64,
641 q: &mut crate::MatrixF64,
642 r: &mut crate::MatrixF64,
643 tau: &mut crate::VectorF64,
644 p: &mut ::Permutation,
645 signum: &mut i32,
646 norm: &mut crate::VectorF64,
647) -> Result<(), Value> {
648 let ret = unsafe {
649 sys::gsl_linalg_QRPT_decomp2(
650 a.unwrap_shared(),
651 q.unwrap_unique(),
652 r.unwrap_unique(),
653 tau.unwrap_unique(),
654 p.unwrap_unique(),
655 signum,
656 norm.unwrap_unique(),
657 )
658 };
659 result_handler!(ret, ())
660}
661
662#[doc(alias = "gsl_linalg_QRPT_solve")]
665pub fn QRPT_solve(
666 qr: &::MatrixF64,
667 tau: &::VectorF64,
668 p: &::Permutation,
669 b: &::VectorF64,
670 x: &mut crate::VectorF64,
671) -> Result<(), Value> {
672 let ret = unsafe {
673 sys::gsl_linalg_QRPT_solve(
674 qr.unwrap_shared(),
675 tau.unwrap_shared(),
676 p.unwrap_shared(),
677 b.unwrap_shared(),
678 x.unwrap_unique(),
679 )
680 };
681 result_handler!(ret, ())
682}
683
684#[doc(alias = "gsl_linalg_QRPT_svx")]
687pub fn QRPT_svx(
688 qr: &::MatrixF64,
689 tau: &::VectorF64,
690 p: &::Permutation,
691 x: &mut crate::VectorF64,
692) -> Result<(), Value> {
693 let ret = unsafe {
694 sys::gsl_linalg_QRPT_svx(
695 qr.unwrap_shared(),
696 tau.unwrap_shared(),
697 p.unwrap_shared(),
698 x.unwrap_unique(),
699 )
700 };
701 result_handler!(ret, ())
702}
703
704#[doc(alias = "gsl_linalg_QRPT_QRsolve")]
707pub fn QRPT_QRsolve(
708 q: &::MatrixF64,
709 r: &::MatrixF64,
710 p: &::Permutation,
711 b: &::VectorF64,
712 x: &mut crate::VectorF64,
713) -> Result<(), Value> {
714 let ret = unsafe {
715 sys::gsl_linalg_QRPT_QRsolve(
716 q.unwrap_shared(),
717 r.unwrap_shared(),
718 p.unwrap_shared(),
719 b.unwrap_shared(),
720 x.unwrap_unique(),
721 )
722 };
723 result_handler!(ret, ())
724}
725
726#[doc(alias = "gsl_linalg_QRPT_update")]
729pub fn QRPT_update(
730 q: &mut crate::MatrixF64,
731 r: &mut crate::MatrixF64,
732 p: &::Permutation,
733 w: &mut crate::VectorF64,
734 v: &::VectorF64,
735) -> Result<(), Value> {
736 let ret = unsafe {
737 sys::gsl_linalg_QRPT_update(
738 q.unwrap_unique(),
739 r.unwrap_unique(),
740 p.unwrap_shared(),
741 w.unwrap_unique(),
742 v.unwrap_shared(),
743 )
744 };
745 result_handler!(ret, ())
746}
747
748#[doc(alias = "gsl_linalg_QRPT_Rsolve")]
750pub fn QRPT_Rsolve(
751 qr: &::MatrixF64,
752 p: &::Permutation,
753 b: &::VectorF64,
754 x: &mut crate::VectorF64,
755) -> Result<(), Value> {
756 let ret = unsafe {
757 sys::gsl_linalg_QRPT_Rsolve(
758 qr.unwrap_shared(),
759 p.unwrap_shared(),
760 b.unwrap_shared(),
761 x.unwrap_unique(),
762 )
763 };
764 result_handler!(ret, ())
765}
766
767#[doc(alias = "gsl_linalg_QRPT_Rsvx")]
770pub fn QRPT_Rsvx(
771 qr: &::MatrixF64,
772 p: &::Permutation,
773 x: &mut crate::VectorF64,
774) -> Result<(), Value> {
775 let ret = unsafe {
776 sys::gsl_linalg_QRPT_Rsvx(qr.unwrap_shared(), p.unwrap_shared(), x.unwrap_unique())
777 };
778 result_handler!(ret, ())
779}
780
781#[doc(alias = "gsl_linalg_SV_decomp")]
788pub fn SV_decomp(
789 a: &mut crate::MatrixF64,
790 v: &mut crate::MatrixF64,
791 s: &mut crate::VectorF64,
792 work: &mut crate::VectorF64,
793) -> Result<(), Value> {
794 let ret = unsafe {
795 sys::gsl_linalg_SV_decomp(
796 a.unwrap_unique(),
797 v.unwrap_unique(),
798 s.unwrap_unique(),
799 work.unwrap_unique(),
800 )
801 };
802 result_handler!(ret, ())
803}
804
805#[doc(alias = "gsl_linalg_SV_decomp_mod")]
808pub fn SV_decomp_mod(
809 a: &mut crate::MatrixF64,
810 x: &mut crate::MatrixF64,
811 v: &mut crate::MatrixF64,
812 s: &mut crate::VectorF64,
813 work: &mut crate::VectorF64,
814) -> Result<(), Value> {
815 let ret = unsafe {
816 sys::gsl_linalg_SV_decomp_mod(
817 a.unwrap_unique(),
818 x.unwrap_unique(),
819 v.unwrap_unique(),
820 s.unwrap_unique(),
821 work.unwrap_unique(),
822 )
823 };
824 result_handler!(ret, ())
825}
826
827#[doc(alias = "gsl_linalg_SV_decomp_jacobi")]
830pub fn SV_decomp_jacobi(
831 a: &mut crate::MatrixF64,
832 v: &mut crate::MatrixF64,
833 s: &mut crate::VectorF64,
834) -> Result<(), Value> {
835 let ret = unsafe {
836 sys::gsl_linalg_SV_decomp_jacobi(a.unwrap_unique(), v.unwrap_unique(), s.unwrap_unique())
837 };
838 result_handler!(ret, ())
839}
840
841#[doc(alias = "gsl_linalg_SV_solve")]
850pub fn SV_solve(
851 u: &::MatrixF64,
852 v: &::MatrixF64,
853 s: &::VectorF64,
854 b: &::VectorF64,
855 x: &mut crate::VectorF64,
856) -> Result<(), Value> {
857 let ret = unsafe {
858 sys::gsl_linalg_SV_solve(
859 u.unwrap_shared(),
860 v.unwrap_shared(),
861 s.unwrap_shared(),
862 b.unwrap_shared(),
863 x.unwrap_unique(),
864 )
865 };
866 result_handler!(ret, ())
867}
868
869#[doc(alias = "gsl_linalg_SV_leverage")]
873pub fn SV_leverage(u: &::MatrixF64, h: &mut crate::VectorF64) -> Result<(), Value> {
874 let ret = unsafe { sys::gsl_linalg_SV_leverage(u.unwrap_shared(), h.unwrap_unique()) };
875 result_handler!(ret, ())
876}
877
878#[doc(alias = "gsl_linalg_cholesky_decomp")]
886pub fn cholesky_decomp(a: &mut crate::MatrixF64) -> Result<(), Value> {
887 let ret = unsafe { sys::gsl_linalg_cholesky_decomp(a.unwrap_unique()) };
888 result_handler!(ret, ())
889}
890
891#[doc(alias = "gsl_linalg_complex_cholesky_decomp")]
899pub fn complex_cholesky_decomp(a: &mut crate::MatrixComplexF64) -> Result<(), Value> {
900 let ret = unsafe { sys::gsl_linalg_complex_cholesky_decomp(a.unwrap_unique()) };
901 result_handler!(ret, ())
902}
903
904#[doc(alias = "gsl_linalg_cholesky_solve")]
907pub fn cholesky_solve(
908 cholesky: &::MatrixF64,
909 b: &::VectorF64,
910 x: &mut crate::VectorF64,
911) -> Result<(), Value> {
912 let ret = unsafe {
913 sys::gsl_linalg_cholesky_solve(
914 cholesky.unwrap_shared(),
915 b.unwrap_shared(),
916 x.unwrap_unique(),
917 )
918 };
919 result_handler!(ret, ())
920}
921
922#[doc(alias = "gsl_linalg_complex_cholesky_solve")]
925pub fn complex_cholesky_solve(
926 cholesky: &::MatrixComplexF64,
927 b: &::VectorComplexF64,
928 x: &mut crate::VectorComplexF64,
929) -> Result<(), Value> {
930 let ret = unsafe {
931 sys::gsl_linalg_complex_cholesky_solve(
932 cholesky.unwrap_shared(),
933 b.unwrap_shared(),
934 x.unwrap_unique(),
935 )
936 };
937 result_handler!(ret, ())
938}
939
940#[doc(alias = "gsl_linalg_cholesky_svx")]
944pub fn cholesky_svx(cholesky: &::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
945 let ret = unsafe { sys::gsl_linalg_cholesky_svx(cholesky.unwrap_shared(), x.unwrap_unique()) };
946 result_handler!(ret, ())
947}
948
949#[doc(alias = "gsl_linalg_complex_cholesky_svx")]
953pub fn complex_cholesky_svx(
954 cholesky: &::MatrixComplexF64,
955 x: &mut crate::VectorComplexF64,
956) -> Result<(), Value> {
957 let ret = unsafe {
958 sys::gsl_linalg_complex_cholesky_svx(cholesky.unwrap_shared(), x.unwrap_unique())
959 };
960 result_handler!(ret, ())
961}
962
963#[doc(alias = "gsl_linalg_cholesky_invert")]
966pub fn cholesky_invert(cholesky: &mut crate::MatrixF64) -> Result<(), Value> {
967 let ret = unsafe { sys::gsl_linalg_cholesky_invert(cholesky.unwrap_unique()) };
968 result_handler!(ret, ())
969}
970
971#[doc(alias = "gsl_linalg_complex_cholesky_invert")]
974pub fn complex_cholesky_invert(cholesky: &mut crate::MatrixComplexF64) -> Result<(), Value> {
975 let ret = unsafe { sys::gsl_linalg_complex_cholesky_invert(cholesky.unwrap_unique()) };
976 result_handler!(ret, ())
977}
978
979#[doc(alias = "gsl_linalg_symmtd_decomp")]
984pub fn symmtd_decomp(a: &mut crate::MatrixF64, tau: &mut crate::VectorF64) -> Result<(), Value> {
985 let ret = unsafe { sys::gsl_linalg_symmtd_decomp(a.unwrap_unique(), tau.unwrap_unique()) };
986 result_handler!(ret, ())
987}
988
989#[doc(alias = "gsl_linalg_symmtd_unpack")]
992pub fn symmtd_unpack(
993 a: &::MatrixF64,
994 tau: &::VectorF64,
995 q: &mut crate::MatrixF64,
996 diag: &mut crate::VectorF64,
997 subdiag: &mut crate::VectorF64,
998) -> Result<(), Value> {
999 let ret = unsafe {
1000 sys::gsl_linalg_symmtd_unpack(
1001 a.unwrap_shared(),
1002 tau.unwrap_shared(),
1003 q.unwrap_unique(),
1004 diag.unwrap_unique(),
1005 subdiag.unwrap_unique(),
1006 )
1007 };
1008 result_handler!(ret, ())
1009}
1010
1011#[doc(alias = "gsl_linalg_symmtd_unpack_T")]
1014pub fn symmtd_unpack_T(
1015 a: &::MatrixF64,
1016 diag: &mut crate::VectorF64,
1017 subdiag: &mut crate::VectorF64,
1018) -> Result<(), Value> {
1019 let ret = unsafe {
1020 sys::gsl_linalg_symmtd_unpack_T(
1021 a.unwrap_shared(),
1022 diag.unwrap_unique(),
1023 subdiag.unwrap_unique(),
1024 )
1025 };
1026 result_handler!(ret, ())
1027}
1028
1029#[doc(alias = "gsl_linalg_hermtd_decomp")]
1034pub fn hermtd_decomp(
1035 a: &mut crate::MatrixComplexF64,
1036 tau: &mut crate::VectorComplexF64,
1037) -> Result<(), Value> {
1038 let ret = unsafe { sys::gsl_linalg_hermtd_decomp(a.unwrap_unique(), tau.unwrap_unique()) };
1039 result_handler!(ret, ())
1040}
1041
1042#[doc(alias = "gsl_linalg_hermtd_unpack")]
1045pub fn hermtd_unpack(
1046 a: &::MatrixComplexF64,
1047 tau: &::VectorComplexF64,
1048 u: &mut crate::MatrixComplexF64,
1049 diag: &mut crate::VectorF64,
1050 subdiag: &mut crate::VectorF64,
1051) -> Result<(), Value> {
1052 let ret = unsafe {
1053 sys::gsl_linalg_hermtd_unpack(
1054 a.unwrap_shared(),
1055 tau.unwrap_shared(),
1056 u.unwrap_unique(),
1057 diag.unwrap_unique(),
1058 subdiag.unwrap_unique(),
1059 )
1060 };
1061 result_handler!(ret, ())
1062}
1063
1064#[doc(alias = "gsl_linalg_hermtd_unpack_T")]
1067pub fn hermtd_unpack_T(
1068 a: &::MatrixComplexF64,
1069 diag: &mut crate::VectorF64,
1070 subdiag: &mut crate::VectorF64,
1071) -> Result<(), Value> {
1072 let ret = unsafe {
1073 sys::gsl_linalg_hermtd_unpack_T(
1074 a.unwrap_shared(),
1075 diag.unwrap_unique(),
1076 subdiag.unwrap_unique(),
1077 )
1078 };
1079 result_handler!(ret, ())
1080}
1081
1082#[doc(alias = "gsl_linalg_hessenberg_decomp")]
1087pub fn hessenberg_decomp(
1088 a: &mut crate::MatrixF64,
1089 tau: &mut crate::VectorF64,
1090) -> Result<(), Value> {
1091 let ret = unsafe { sys::gsl_linalg_hessenberg_decomp(a.unwrap_unique(), tau.unwrap_unique()) };
1092 result_handler!(ret, ())
1093}
1094
1095#[doc(alias = "gsl_linalg_hessenberg_unpack")]
1098pub fn hessenberg_unpack(
1099 h: &mut crate::MatrixF64,
1100 tau: &mut crate::VectorF64,
1101 u: &mut crate::MatrixF64,
1102) -> Result<(), Value> {
1103 let ret = unsafe {
1104 sys::gsl_linalg_hessenberg_unpack(h.unwrap_unique(), tau.unwrap_unique(), u.unwrap_unique())
1105 };
1106 result_handler!(ret, ())
1107}
1108
1109#[doc(alias = "gsl_linalg_hessenberg_unpack_accum")]
1113pub fn hessenberg_unpack_accum(
1114 h: &mut crate::MatrixF64,
1115 tau: &mut crate::VectorF64,
1116 v: &mut crate::MatrixF64,
1117) -> Result<(), Value> {
1118 let ret = unsafe {
1119 sys::gsl_linalg_hessenberg_unpack_accum(
1120 h.unwrap_unique(),
1121 tau.unwrap_unique(),
1122 v.unwrap_unique(),
1123 )
1124 };
1125 result_handler!(ret, ())
1126}
1127
1128#[doc(alias = "gsl_linalg_hessenberg_set_zero")]
1131pub fn hessenberg_set_zero(h: &mut crate::MatrixF64) -> Result<(), Value> {
1132 let ret = unsafe { sys::gsl_linalg_hessenberg_set_zero(h.unwrap_unique()) };
1133 result_handler!(ret, ())
1134}
1135
1136#[doc(alias = "gsl_linalg_hesstri_decomp")]
1140pub fn hesstri_decomp(
1141 a: &mut crate::MatrixF64,
1142 b: &mut crate::MatrixF64,
1143 u: &mut crate::MatrixF64,
1144 v: &mut crate::MatrixF64,
1145 work: &mut crate::VectorF64,
1146) -> Result<(), Value> {
1147 let ret = unsafe {
1148 sys::gsl_linalg_hesstri_decomp(
1149 a.unwrap_unique(),
1150 b.unwrap_unique(),
1151 u.unwrap_unique(),
1152 v.unwrap_unique(),
1153 work.unwrap_unique(),
1154 )
1155 };
1156 result_handler!(ret, ())
1157}
1158
1159#[doc(alias = "gsl_linalg_bidiag_decomp")]
1164pub fn bidiag_decomp(
1165 a: &mut crate::MatrixF64,
1166 tau_u: &mut crate::VectorF64,
1167 tau_v: &mut crate::VectorF64,
1168) -> Result<(), Value> {
1169 let ret = unsafe {
1170 sys::gsl_linalg_bidiag_decomp(
1171 a.unwrap_unique(),
1172 tau_u.unwrap_unique(),
1173 tau_v.unwrap_unique(),
1174 )
1175 };
1176 result_handler!(ret, ())
1177}
1178
1179#[doc(alias = "gsl_linalg_bidiag_unpack")]
1183pub fn bidiag_unpack(
1184 a: &mut crate::MatrixF64,
1185 tau_u: &::VectorF64,
1186 u: &mut crate::MatrixF64,
1187 tau_v: &::VectorF64,
1188 v: &mut crate::MatrixF64,
1189 diag: &mut crate::VectorF64,
1190 superdiag: &mut crate::VectorF64,
1191) -> Result<(), Value> {
1192 let ret = unsafe {
1193 sys::gsl_linalg_bidiag_unpack(
1194 a.unwrap_unique(),
1195 tau_u.unwrap_shared(),
1196 u.unwrap_unique(),
1197 tau_v.unwrap_shared(),
1198 v.unwrap_unique(),
1199 diag.unwrap_unique(),
1200 superdiag.unwrap_unique(),
1201 )
1202 };
1203 result_handler!(ret, ())
1204}
1205
1206#[doc(alias = "gsl_linalg_bidiag_unpack2")]
1209pub fn bidiag_unpack2(
1210 a: &mut crate::MatrixF64,
1211 tau_u: &mut crate::VectorF64,
1212 tau_v: &mut crate::VectorF64,
1213 v: &mut crate::MatrixF64,
1214) -> Result<(), Value> {
1215 let ret = unsafe {
1216 sys::gsl_linalg_bidiag_unpack2(
1217 a.unwrap_unique(),
1218 tau_u.unwrap_unique(),
1219 tau_v.unwrap_unique(),
1220 v.unwrap_unique(),
1221 )
1222 };
1223 result_handler!(ret, ())
1224}
1225
1226#[doc(alias = "gsl_linalg_bidiag_unpack_B")]
1229pub fn bidiag_unpack_B(
1230 a: &::MatrixF64,
1231 diag: &mut crate::VectorF64,
1232 superdiag: &mut crate::VectorF64,
1233) -> Result<(), Value> {
1234 let ret = unsafe {
1235 sys::gsl_linalg_bidiag_unpack_B(
1236 a.unwrap_shared(),
1237 diag.unwrap_unique(),
1238 superdiag.unwrap_unique(),
1239 )
1240 };
1241 result_handler!(ret, ())
1242}
1243
1244#[doc(alias = "gsl_linalg_householder_transform")]
1247pub fn householder_transform(v: &mut crate::VectorF64) -> f64 {
1248 unsafe { sys::gsl_linalg_householder_transform(v.unwrap_unique()) }
1249}
1250
1251#[doc(alias = "gsl_linalg_complex_householder_transform")]
1254pub fn complex_householder_transform(v: &mut crate::VectorComplexF64) -> crate::ComplexF64 {
1255 unsafe {
1256 std::mem::transmute(sys::gsl_linalg_complex_householder_transform(
1257 v.unwrap_unique(),
1258 ))
1259 }
1260}
1261
1262#[doc(alias = "gsl_linalg_householder_hm")]
1265pub fn householder_hm(tau: f64, v: &::VectorF64, a: &mut crate::MatrixF64) -> Result<(), Value> {
1266 let ret = unsafe { sys::gsl_linalg_householder_hm(tau, v.unwrap_shared(), a.unwrap_unique()) };
1267 result_handler!(ret, ())
1268}
1269
1270#[doc(alias = "gsl_linalg_complex_householder_hm")]
1273pub fn complex_householder_hm(
1274 tau: &::ComplexF64,
1275 v: &::VectorComplexF64,
1276 a: &mut crate::MatrixComplexF64,
1277) -> Result<(), Value> {
1278 let ret = unsafe {
1279 sys::gsl_linalg_complex_householder_hm(
1280 std::mem::transmute(*tau),
1281 v.unwrap_shared(),
1282 a.unwrap_unique(),
1283 )
1284 };
1285 result_handler!(ret, ())
1286}
1287
1288#[doc(alias = "gsl_linalg_householder_mh")]
1291pub fn householder_mh(tau: f64, v: &::VectorF64, a: &mut crate::MatrixF64) -> Result<(), Value> {
1292 let ret = unsafe { sys::gsl_linalg_householder_mh(tau, v.unwrap_shared(), a.unwrap_unique()) };
1293 result_handler!(ret, ())
1294}
1295
1296#[doc(alias = "gsl_linalg_complex_householder_mh")]
1299pub fn complex_householder_mh(
1300 tau: &::ComplexF64,
1301 v: &::VectorComplexF64,
1302 a: &mut crate::MatrixComplexF64,
1303) -> Result<(), Value> {
1304 let ret = unsafe {
1305 sys::gsl_linalg_complex_householder_mh(
1306 std::mem::transmute(*tau),
1307 v.unwrap_shared(),
1308 a.unwrap_unique(),
1309 )
1310 };
1311 result_handler!(ret, ())
1312}
1313
1314#[doc(alias = "gsl_linalg_householder_hv")]
1317pub fn householder_hv(tau: f64, v: &::VectorF64, w: &mut crate::VectorF64) -> Result<(), Value> {
1318 let ret = unsafe { sys::gsl_linalg_householder_hv(tau, v.unwrap_shared(), w.unwrap_unique()) };
1319 result_handler!(ret, ())
1320}
1321
1322#[doc(alias = "gsl_linalg_complex_householder_hv")]
1325pub fn complex_householder_hv(
1326 tau: &::ComplexF64,
1327 v: &::VectorComplexF64,
1328 w: &mut crate::VectorComplexF64,
1329) -> Result<(), Value> {
1330 let ret = unsafe {
1331 sys::gsl_linalg_complex_householder_hv(
1332 std::mem::transmute(*tau),
1333 v.unwrap_shared(),
1334 w.unwrap_unique(),
1335 )
1336 };
1337 result_handler!(ret, ())
1338}
1339
1340#[doc(alias = "gsl_linalg_HH_solve")]
1343pub fn HH_solve(
1344 mut a: crate::MatrixF64,
1345 b: &::VectorF64,
1346 x: &mut crate::VectorF64,
1347) -> Result<(), Value> {
1348 let ret = unsafe {
1349 sys::gsl_linalg_HH_solve(a.unwrap_unique(), b.unwrap_shared(), x.unwrap_unique())
1350 };
1351 result_handler!(ret, ())
1352}
1353
1354#[doc(alias = "gsl_linalg_HH_svx")]
1357pub fn HH_svx(mut a: crate::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
1358 let ret = unsafe { sys::gsl_linalg_HH_svx(a.unwrap_unique(), x.unwrap_unique()) };
1359 result_handler!(ret, ())
1360}
1361
1362#[doc(alias = "gsl_linalg_solve_tridiag")]
1372pub fn solve_tridiag(
1373 diag: &::VectorF64,
1374 e: &::VectorF64,
1375 f: &::VectorF64,
1376 b: &::VectorF64,
1377 x: &mut crate::VectorF64,
1378) -> Result<(), Value> {
1379 let ret = unsafe {
1380 sys::gsl_linalg_solve_tridiag(
1381 diag.unwrap_shared(),
1382 e.unwrap_shared(),
1383 f.unwrap_shared(),
1384 b.unwrap_shared(),
1385 x.unwrap_unique(),
1386 )
1387 };
1388 result_handler!(ret, ())
1389}
1390
1391#[doc(alias = "gsl_linalg_solve_symm_tridiag")]
1401pub fn solve_symm_tridiag(
1402 diag: &::VectorF64,
1403 e: &::VectorF64,
1404 b: &::VectorF64,
1405 x: &mut crate::VectorF64,
1406) -> Result<(), Value> {
1407 let ret = unsafe {
1408 sys::gsl_linalg_solve_symm_tridiag(
1409 diag.unwrap_shared(),
1410 e.unwrap_shared(),
1411 b.unwrap_shared(),
1412 x.unwrap_unique(),
1413 )
1414 };
1415 result_handler!(ret, ())
1416}
1417
1418#[doc(alias = "gsl_linalg_solve_cyc_tridiag")]
1428pub fn solve_cyc_tridiag(
1429 diag: &::VectorF64,
1430 e: &::VectorF64,
1431 f: &::VectorF64,
1432 b: &::VectorF64,
1433 x: &mut crate::VectorF64,
1434) -> Result<(), Value> {
1435 let ret = unsafe {
1436 sys::gsl_linalg_solve_cyc_tridiag(
1437 diag.unwrap_shared(),
1438 e.unwrap_shared(),
1439 f.unwrap_shared(),
1440 b.unwrap_shared(),
1441 x.unwrap_unique(),
1442 )
1443 };
1444 result_handler!(ret, ())
1445}
1446
1447#[doc(alias = "gsl_linalg_solve_symm_cyc_tridiag")]
1457pub fn solve_symm_cyc_tridiag(
1458 diag: &::VectorF64,
1459 e: &::VectorF64,
1460 b: &::VectorF64,
1461 x: &mut crate::VectorF64,
1462) -> Result<(), Value> {
1463 let ret = unsafe {
1464 sys::gsl_linalg_solve_symm_cyc_tridiag(
1465 diag.unwrap_shared(),
1466 e.unwrap_shared(),
1467 b.unwrap_shared(),
1468 x.unwrap_unique(),
1469 )
1470 };
1471 result_handler!(ret, ())
1472}
1473
1474#[doc(alias = "gsl_linalg_balance_matrix")]
1477pub fn balance_matrix(a: &mut crate::MatrixF64, d: &mut crate::VectorF64) -> Result<(), Value> {
1478 let ret = unsafe { sys::gsl_linalg_balance_matrix(a.unwrap_unique(), d.unwrap_unique()) };
1479 result_handler!(ret, ())
1480}
1481
1482#[cfg(feature = "v2_2")]
1483#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1484#[doc(alias = "gsl_linalg_pcholesky_decomp")]
1485pub fn pcholesky_decomp(a: &mut crate::MatrixF64, p: &mut ::Permutation) -> Result<(), Value> {
1486 let ret = unsafe { sys::gsl_linalg_pcholesky_decomp(a.unwrap_unique(), p.unwrap_unique()) };
1487 result_handler!(ret, ())
1488}
1489
1490#[cfg(feature = "v2_2")]
1491#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1492#[doc(alias = "gsl_linalg_pcholesky_solve")]
1493pub fn pcholesky_solve(
1494 LDLT: &::MatrixF64,
1495 p: &::Permutation,
1496 b: &::VectorF64,
1497 x: &mut crate::VectorF64,
1498) -> Result<(), Value> {
1499 let ret = unsafe {
1500 sys::gsl_linalg_pcholesky_solve(
1501 LDLT.unwrap_shared(),
1502 p.unwrap_shared(),
1503 b.unwrap_shared(),
1504 x.unwrap_unique(),
1505 )
1506 };
1507 result_handler!(ret, ())
1508}
1509
1510#[cfg(feature = "v2_2")]
1511#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1512#[doc(alias = "gsl_linalg_pcholesky_svx")]
1513pub fn pcholesky_svx(
1514 LDLT: &::MatrixF64,
1515 p: &::Permutation,
1516 x: &mut crate::VectorF64,
1517) -> Result<(), Value> {
1518 let ret = unsafe {
1519 sys::gsl_linalg_pcholesky_svx(LDLT.unwrap_shared(), p.unwrap_shared(), x.unwrap_unique())
1520 };
1521 result_handler!(ret, ())
1522}
1523
1524#[cfg(feature = "v2_2")]
1525#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1526#[doc(alias = "gsl_linalg_pcholesky_decomp2")]
1527pub fn pcholesky_decomp2(
1528 A: &mut crate::MatrixF64,
1529 p: &mut ::Permutation,
1530 S: &mut crate::VectorF64,
1531) -> Result<(), Value> {
1532 let ret = unsafe {
1533 sys::gsl_linalg_pcholesky_decomp2(A.unwrap_unique(), p.unwrap_unique(), S.unwrap_unique())
1534 };
1535 result_handler!(ret, ())
1536}
1537
1538#[cfg(feature = "v2_2")]
1539#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1540#[doc(alias = "gsl_linalg_pcholesky_solve2")]
1541pub fn pcholesky_solve2(
1542 LDLT: &::MatrixF64,
1543 p: &::Permutation,
1544 S: &::VectorF64,
1545 b: &::VectorF64,
1546 x: &mut crate::VectorF64,
1547) -> Result<(), Value> {
1548 let ret = unsafe {
1549 sys::gsl_linalg_pcholesky_solve2(
1550 LDLT.unwrap_shared(),
1551 p.unwrap_shared(),
1552 S.unwrap_shared(),
1553 b.unwrap_shared(),
1554 x.unwrap_unique(),
1555 )
1556 };
1557 result_handler!(ret, ())
1558}
1559
1560#[cfg(feature = "v2_2")]
1561#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1562#[doc(alias = "gsl_linalg_pcholesky_svx2")]
1563pub fn pcholesky_svx2(
1564 LDLT: &::MatrixF64,
1565 p: &::Permutation,
1566 S: &::VectorF64,
1567 x: &mut crate::VectorF64,
1568) -> Result<(), Value> {
1569 let ret = unsafe {
1570 sys::gsl_linalg_pcholesky_svx2(
1571 LDLT.unwrap_shared(),
1572 p.unwrap_shared(),
1573 S.unwrap_shared(),
1574 x.unwrap_unique(),
1575 )
1576 };
1577 result_handler!(ret, ())
1578}
1579
1580#[cfg(feature = "v2_2")]
1581#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1582#[doc(alias = "gsl_linalg_pcholesky_invert")]
1583pub fn pcholesky_invert(
1584 LDLT: &::MatrixF64,
1585 p: &::Permutation,
1586 Ainv: &mut crate::MatrixF64,
1587) -> Result<(), Value> {
1588 let ret = unsafe {
1589 sys::gsl_linalg_pcholesky_invert(
1590 LDLT.unwrap_shared(),
1591 p.unwrap_shared(),
1592 Ainv.unwrap_unique(),
1593 )
1594 };
1595 result_handler!(ret, ())
1596}
1597
1598#[cfg(feature = "v2_2")]
1600#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1601#[doc(alias = "gsl_linalg_pcholesky_rcond")]
1602pub fn pcholesky_rcond(
1603 LDLT: &::MatrixF64,
1604 p: &::Permutation,
1605 work: &mut crate::VectorF64,
1606) -> Result<f64, Value> {
1607 let mut rcond = 0.;
1608 let ret = unsafe {
1609 sys::gsl_linalg_pcholesky_rcond(
1610 LDLT.unwrap_shared(),
1611 p.unwrap_shared(),
1612 &mut rcond,
1613 work.unwrap_unique(),
1614 )
1615 };
1616 result_handler!(ret, rcond)
1617}
1618
1619#[cfg(feature = "v2_2")]
1620#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1621#[doc(alias = "gsl_linalg_mcholesky_decomp")]
1622pub fn mcholesky_decomp(
1623 A: &mut crate::MatrixF64,
1624 p: &mut ::Permutation,
1625 E: &mut crate::VectorF64,
1626) -> Result<(), Value> {
1627 let ret = unsafe {
1628 sys::gsl_linalg_mcholesky_decomp(A.unwrap_unique(), p.unwrap_unique(), E.unwrap_unique())
1629 };
1630 result_handler!(ret, ())
1631}
1632
1633#[cfg(feature = "v2_2")]
1634#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1635#[doc(alias = "gsl_linalg_mcholesky_solve")]
1636pub fn mcholesky_solve(
1637 LDLT: &::MatrixF64,
1638 p: &::Permutation,
1639 b: &::VectorF64,
1640 x: &mut crate::VectorF64,
1641) -> Result<(), Value> {
1642 let ret = unsafe {
1643 sys::gsl_linalg_mcholesky_solve(
1644 LDLT.unwrap_shared(),
1645 p.unwrap_shared(),
1646 b.unwrap_shared(),
1647 x.unwrap_unique(),
1648 )
1649 };
1650 result_handler!(ret, ())
1651}
1652
1653#[cfg(feature = "v2_2")]
1654#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1655#[doc(alias = "gsl_linalg_mcholesky_svx")]
1656pub fn mcholesky_svx(
1657 LDLT: &::MatrixF64,
1658 p: &::Permutation,
1659 x: &mut crate::VectorF64,
1660) -> Result<(), Value> {
1661 let ret = unsafe {
1662 sys::gsl_linalg_mcholesky_svx(LDLT.unwrap_shared(), p.unwrap_shared(), x.unwrap_unique())
1663 };
1664 result_handler!(ret, ())
1665}
1666
1667#[cfg(feature = "v2_2")]
1669#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1670#[doc(alias = "gsl_linalg_mcholesky_rcond")]
1671pub fn mcholesky_rcond(
1672 LDLT: &::MatrixF64,
1673 p: &::Permutation,
1674 work: &mut crate::VectorF64,
1675) -> Result<f64, Value> {
1676 let mut rcond = 0.;
1677 let ret = unsafe {
1678 sys::gsl_linalg_mcholesky_rcond(
1679 LDLT.unwrap_shared(),
1680 p.unwrap_shared(),
1681 &mut rcond,
1682 work.unwrap_unique(),
1683 )
1684 };
1685 result_handler!(ret, rcond)
1686}
1687
1688#[cfg(feature = "v2_2")]
1689#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1690#[doc(alias = "gsl_linalg_mcholesky_invert")]
1691pub fn mcholesky_invert(
1692 LDLT: &::MatrixF64,
1693 p: &::Permutation,
1694 Ainv: &mut crate::MatrixF64,
1695) -> Result<(), Value> {
1696 let ret = unsafe {
1697 sys::gsl_linalg_mcholesky_invert(
1698 LDLT.unwrap_shared(),
1699 p.unwrap_shared(),
1700 Ainv.unwrap_unique(),
1701 )
1702 };
1703 result_handler!(ret, ())
1704}
1705
1706#[cfg(feature = "v2_6")]
1707#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1708#[doc(alias = "gsl_linalg_cholesky_band_decomp")]
1709pub fn cholesky_band_decomp(A: &mut crate::MatrixF64) -> Result<(), Value> {
1710 let ret = unsafe { sys::gsl_linalg_cholesky_band_decomp(A.unwrap_unique()) };
1711 result_handler!(ret, ())
1712}
1713
1714#[cfg(feature = "v2_6")]
1715#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1716#[doc(alias = "gsl_linalg_cholesky_band_solve")]
1717pub fn cholesky_band_solve(
1718 LLT: &::MatrixF64,
1719 b: &::VectorF64,
1720 x: &mut crate::VectorF64,
1721) -> Result<(), Value> {
1722 let ret = unsafe {
1723 sys::gsl_linalg_cholesky_band_solve(
1724 LLT.unwrap_shared(),
1725 b.unwrap_shared(),
1726 x.unwrap_unique(),
1727 )
1728 };
1729 result_handler!(ret, ())
1730}
1731
1732#[cfg(feature = "v2_6")]
1733#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1734#[doc(alias = "gsl_linalg_cholesky_band_svx")]
1735pub fn cholesky_band_svx(LLT: &::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
1736 let ret = unsafe { sys::gsl_linalg_cholesky_band_svx(LLT.unwrap_shared(), x.unwrap_unique()) };
1737 result_handler!(ret, ())
1738}
1739
1740#[cfg(feature = "v2_7")]
1741#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_7")))]
1742#[doc(alias = "gsl_linalg_cholesky_band_solvem")]
1743pub fn cholesky_band_solvem(
1744 LLT: &::MatrixF64,
1745 B: &::MatrixF64,
1746 X: &mut crate::MatrixF64,
1747) -> Result<(), Value> {
1748 let ret = unsafe {
1749 sys::gsl_linalg_cholesky_band_solvem(
1750 LLT.unwrap_shared(),
1751 B.unwrap_shared(),
1752 X.unwrap_unique(),
1753 )
1754 };
1755 result_handler!(ret, ())
1756}
1757
1758#[cfg(feature = "v2_7")]
1759#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_7")))]
1760#[doc(alias = "gsl_linalg_cholesky_band_svxm")]
1761pub fn cholesky_band_svxm(LLT: &::MatrixF64, X: &mut crate::MatrixF64) -> Result<(), Value> {
1762 let ret = unsafe { sys::gsl_linalg_cholesky_band_svxm(LLT.unwrap_shared(), X.unwrap_unique()) };
1763 result_handler!(ret, ())
1764}
1765
1766#[cfg(feature = "v2_6")]
1767#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1768#[doc(alias = "gsl_linalg_cholesky_band_invert")]
1769pub fn cholesky_band_invert(LLT: &::MatrixF64, Ainv: &mut crate::MatrixF64) -> Result<(), Value> {
1770 let ret =
1771 unsafe { sys::gsl_linalg_cholesky_band_invert(LLT.unwrap_shared(), Ainv.unwrap_unique()) };
1772 result_handler!(ret, ())
1773}
1774
1775#[cfg(feature = "v2_6")]
1776#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1777#[doc(alias = "gsl_linalg_cholesky_band_unpack")]
1778pub fn cholesky_band_unpack(LLT: &::MatrixF64, L: &mut crate::MatrixF64) -> Result<(), Value> {
1779 let ret =
1780 unsafe { sys::gsl_linalg_cholesky_band_unpack(LLT.unwrap_shared(), L.unwrap_unique()) };
1781 result_handler!(ret, ())
1782}
1783
1784#[cfg(feature = "v2_6")]
1786#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1787#[doc(alias = "gsl_linalg_cholesky_band_rcond")]
1788pub fn cholesky_band_rcond(LLT: &::MatrixF64, work: &mut crate::VectorF64) -> Result<f64, Value> {
1789 let mut rcond = 0.;
1790 let ret = unsafe {
1791 sys::gsl_linalg_cholesky_band_rcond(LLT.unwrap_shared(), &mut rcond, work.unwrap_unique())
1792 };
1793 result_handler!(ret, rcond)
1794}
1795
1796#[cfg(feature = "v2_6")]
1797#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1798#[doc(alias = "gsl_linalg_ldlt_decomp")]
1799pub fn ldlt_decomp(A: &mut crate::MatrixF64) -> Result<(), Value> {
1800 let ret = unsafe { sys::gsl_linalg_ldlt_decomp(A.unwrap_unique()) };
1801 result_handler!(ret, ())
1802}
1803
1804#[cfg(feature = "v2_6")]
1805#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1806#[doc(alias = "gsl_linalg_ldlt_solve")]
1807pub fn ldlt_solve(
1808 LDLT: &::MatrixF64,
1809 b: &::VectorF64,
1810 x: &mut crate::VectorF64,
1811) -> Result<(), Value> {
1812 let ret = unsafe {
1813 sys::gsl_linalg_ldlt_solve(LDLT.unwrap_shared(), b.unwrap_shared(), x.unwrap_unique())
1814 };
1815 result_handler!(ret, ())
1816}
1817
1818#[cfg(feature = "v2_6")]
1819#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1820#[doc(alias = "gsl_linalg_ldlt_svx")]
1821pub fn ldlt_svx(LDLT: &::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
1822 let ret = unsafe { sys::gsl_linalg_ldlt_svx(LDLT.unwrap_shared(), x.unwrap_unique()) };
1823 result_handler!(ret, ())
1824}
1825
1826#[cfg(feature = "v2_6")]
1828#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1829#[doc(alias = "gsl_linalg_ldlt_rcond")]
1830pub fn ldlt_rcond(LDLT: &::MatrixF64, work: &mut crate::VectorF64) -> Result<f64, Value> {
1831 let mut rcond = 0.;
1832 let ret = unsafe {
1833 sys::gsl_linalg_ldlt_rcond(LDLT.unwrap_shared(), &mut rcond, work.unwrap_unique())
1834 };
1835 result_handler!(ret, rcond)
1836}
1837
1838#[cfg(feature = "v2_6")]
1839#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1840#[doc(alias = "gsl_linalg_ldlt_band_decomp")]
1841pub fn ldlt_band_decomp(A: &mut crate::MatrixF64) -> Result<(), Value> {
1842 let ret = unsafe { sys::gsl_linalg_ldlt_band_decomp(A.unwrap_unique()) };
1843 result_handler!(ret, ())
1844}
1845
1846#[cfg(feature = "v2_6")]
1847#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1848#[doc(alias = "gsl_linalg_ldlt_band_solve")]
1849pub fn ldlt_band_solve(
1850 LDLT: &::MatrixF64,
1851 b: &::VectorF64,
1852 x: &mut crate::VectorF64,
1853) -> Result<(), Value> {
1854 let ret = unsafe {
1855 sys::gsl_linalg_ldlt_band_solve(LDLT.unwrap_shared(), b.unwrap_shared(), x.unwrap_unique())
1856 };
1857 result_handler!(ret, ())
1858}
1859
1860#[cfg(feature = "v2_6")]
1861#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1862#[doc(alias = "gsl_linalg_ldlt_band_svx")]
1863pub fn ldlt_band_svx(LDLT: &::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
1864 let ret = unsafe { sys::gsl_linalg_ldlt_band_svx(LDLT.unwrap_shared(), x.unwrap_unique()) };
1865 result_handler!(ret, ())
1866}
1867
1868#[cfg(feature = "v2_6")]
1869#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1870#[doc(alias = "gsl_linalg_ldlt_band_unpack")]
1871pub fn ldlt_band_unpack(
1872 LDLT: &::MatrixF64,
1873 L: &mut crate::MatrixF64,
1874 D: &mut crate::VectorF64,
1875) -> Result<(), Value> {
1876 let ret = unsafe {
1877 sys::gsl_linalg_ldlt_band_unpack(LDLT.unwrap_shared(), L.unwrap_unique(), D.unwrap_unique())
1878 };
1879 result_handler!(ret, ())
1880}
1881
1882#[cfg(feature = "v2_6")]
1884#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1885#[doc(alias = "gsl_linalg_ldlt_band_rcond")]
1886pub fn ldlt_band_rcond(LDLT: &::MatrixF64, work: &mut crate::VectorF64) -> Result<f64, Value> {
1887 let mut rcond = 0.;
1888 let ret = unsafe {
1889 sys::gsl_linalg_ldlt_band_rcond(LDLT.unwrap_shared(), &mut rcond, work.unwrap_unique())
1890 };
1891 result_handler!(ret, rcond)
1892}
1893
1894#[cfg(feature = "v2_2")]
1895#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1896#[doc(alias = "gsl_linalg_tri_upper_invert")]
1897pub fn tri_upper_invert(T: &mut crate::MatrixF64) -> Result<(), Value> {
1898 let ret = unsafe { sys::gsl_linalg_tri_upper_invert(T.unwrap_unique()) };
1899 result_handler!(ret, ())
1900}
1901
1902#[cfg(feature = "v2_2")]
1903#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1904#[doc(alias = "gsl_linalg_tri_lower_invert")]
1905pub fn tri_lower_invert(T: &mut crate::MatrixF64) -> Result<(), Value> {
1906 let ret = unsafe { sys::gsl_linalg_tri_lower_invert(T.unwrap_unique()) };
1907 result_handler!(ret, ())
1908}
1909
1910#[cfg(feature = "v2_2")]
1911#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1912#[doc(alias = "gsl_linalg_tri_upper_unit_invert")]
1913pub fn tri_upper_unit_invert(T: &mut crate::MatrixF64) -> Result<(), Value> {
1914 let ret = unsafe { sys::gsl_linalg_tri_upper_unit_invert(T.unwrap_unique()) };
1915 result_handler!(ret, ())
1916}
1917
1918#[cfg(feature = "v2_2")]
1919#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1920#[doc(alias = "gsl_linalg_tri_lower_unit_invert")]
1921pub fn tri_lower_unit_invert(T: &mut crate::MatrixF64) -> Result<(), Value> {
1922 let ret = unsafe { sys::gsl_linalg_tri_lower_unit_invert(T.unwrap_unique()) };
1923 result_handler!(ret, ())
1924}
1925
1926#[cfg(feature = "v2_2")]
1927#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1928#[doc(alias = "gsl_linalg_complex_tri_invert")]
1929pub fn tri_invert(
1930 Uplo: enums::CblasUplo,
1931 Diag: enums::CblasDiag,
1932 T: &mut crate::MatrixComplexF64,
1933) -> Result<(), Value> {
1934 let ret =
1935 unsafe { sys::gsl_linalg_complex_tri_invert(Uplo.into(), Diag.into(), T.unwrap_unique()) };
1936 result_handler!(ret, ())
1937}
1938
1939#[doc(alias = "gsl_linalg_complex_tri_invert")]
1940pub fn complex_tri_invert(
1941 Uplo: enums::CblasUplo,
1942 Diag: enums::CblasDiag,
1943 T: &mut crate::MatrixComplexF64,
1944) -> Result<(), Value> {
1945 let ret =
1946 unsafe { sys::gsl_linalg_complex_tri_invert(Uplo.into(), Diag.into(), T.unwrap_unique()) };
1947 result_handler!(ret, ())
1948}
1949
1950#[cfg(feature = "v2_2")]
1951#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1952#[doc(alias = "gsl_linalg_tri_LTL")]
1953pub fn tri_LTL(L: &mut crate::MatrixF64) -> Result<(), Value> {
1954 let ret = unsafe { sys::gsl_linalg_tri_LTL(L.unwrap_unique()) };
1955 result_handler!(ret, ())
1956}
1957
1958#[cfg(feature = "v2_2")]
1959#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1960#[doc(alias = "gsl_linalg_tri_UL")]
1961pub fn tri_UL(LU: &mut crate::MatrixF64) -> Result<(), Value> {
1962 let ret = unsafe { sys::gsl_linalg_tri_UL(LU.unwrap_unique()) };
1963 result_handler!(ret, ())
1964}
1965
1966#[doc(alias = "gsl_linalg_complex_tri_LHL")]
1967pub fn complex_tri_LHL(L: &mut crate::MatrixComplexF64) -> Result<(), Value> {
1968 let ret = unsafe { sys::gsl_linalg_complex_tri_LHL(L.unwrap_unique()) };
1969 result_handler!(ret, ())
1970}
1971
1972#[doc(alias = "gsl_linalg_complex_tri_UL")]
1973pub fn complex_tri_UL(LU: &mut crate::MatrixComplexF64) -> Result<(), Value> {
1974 let ret = unsafe { sys::gsl_linalg_complex_tri_UL(LU.unwrap_unique()) };
1975 result_handler!(ret, ())
1976}
1977
1978#[doc(alias = "gsl_linalg_givens")]
1980pub fn givens(a: f64, b: f64) -> (f64, f64) {
1981 let mut c = 0.;
1982 let mut s = 0.;
1983 unsafe { sys::gsl_linalg_givens(a, b, &mut c, &mut s) };
1984 (c, s)
1985}
1986
1987#[doc(alias = "gsl_linalg_givens_gv")]
1988pub fn givens_gv(v: &mut crate::VectorF64, i: usize, j: usize, c: f64, s: f64) {
1989 unsafe { sys::gsl_linalg_givens_gv(v.unwrap_unique(), i, j, c, s) }
1990}