Free cookie consent management tool by TermsFeed Policy Generator

source: branches/3.2/sources/HeuristicLab.ExtLibs/HeuristicLab.ALGLIB/2.1.2.2591/ALGLIB-2.1.2.2591/sdet.cs @ 4539

Last change on this file since 4539 was 2645, checked in by mkommend, 15 years ago

extracted external libraries and adapted dependent plugins (ticket #837)

File size: 5.8 KB
Line 
1/*************************************************************************
2Copyright (c) 2005-2007, Sergey Bochkanov (ALGLIB project).
3
4>>> SOURCE LICENSE >>>
5This program is free software; you can redistribute it and/or modify
6it under the terms of the GNU General Public License as published by
7the Free Software Foundation (www.fsf.org); either version 2 of the
8License, or (at your option) any later version.
9
10This program is distributed in the hope that it will be useful,
11but WITHOUT ANY WARRANTY; without even the implied warranty of
12MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
13GNU General Public License for more details.
14
15A copy of the GNU General Public License is available at
16http://www.fsf.org/licensing/licenses
17
18>>> END OF LICENSE >>>
19*************************************************************************/
20
21using System;
22
23namespace alglib
24{
25    public class sdet
26    {
27        /*************************************************************************
28        Determinant calculation of the matrix given by LDLT decomposition.
29
30        Input parameters:
31            A       -   LDLT-decomposition of the matrix,
32                        output of subroutine SMatrixLDLT.
33            Pivots  -   table of permutations which were made during
34                        LDLT decomposition, output of subroutine SMatrixLDLT.
35            N       -   size of matrix A.
36            IsUpper -   matrix storage format. The value is equal to the input
37                        parameter of subroutine SMatrixLDLT.
38
39        Result:
40            matrix determinant.
41
42          -- ALGLIB --
43             Copyright 2005-2008 by Bochkanov Sergey
44        *************************************************************************/
45        public static double smatrixldltdet(ref double[,] a,
46            ref int[] pivots,
47            int n,
48            bool isupper)
49        {
50            double result = 0;
51            int k = 0;
52
53            result = 1;
54            if( isupper )
55            {
56                k = 0;
57                while( k<n )
58                {
59                    if( pivots[k]>=0 )
60                    {
61                        result = result*a[k,k];
62                        k = k+1;
63                    }
64                    else
65                    {
66                        result = result*(a[k,k]*a[k+1,k+1]-a[k,k+1]*a[k,k+1]);
67                        k = k+2;
68                    }
69                }
70            }
71            else
72            {
73                k = n-1;
74                while( k>=0 )
75                {
76                    if( pivots[k]>=0 )
77                    {
78                        result = result*a[k,k];
79                        k = k-1;
80                    }
81                    else
82                    {
83                        result = result*(a[k-1,k-1]*a[k,k]-a[k,k-1]*a[k,k-1]);
84                        k = k-2;
85                    }
86                }
87            }
88            return result;
89        }
90
91
92        /*************************************************************************
93        Determinant calculation of the symmetric matrix
94
95        Input parameters:
96            A       -   matrix. Array with elements [0..N-1, 0..N-1].
97            N       -   size of matrix A.
98            IsUpper -   if IsUpper = True, then symmetric matrix A is given by its
99                        upper triangle, and the lower triangle isn’t used by
100                        subroutine. Similarly, if IsUpper = False, then A is given
101                        by its lower triangle.
102
103        Result:
104            determinant of matrix A.
105
106          -- ALGLIB --
107             Copyright 2005-2008 by Bochkanov Sergey
108        *************************************************************************/
109        public static double smatrixdet(double[,] a,
110            int n,
111            bool isupper)
112        {
113            double result = 0;
114            int[] pivots = new int[0];
115
116            a = (double[,])a.Clone();
117
118            ldlt.smatrixldlt(ref a, n, isupper, ref pivots);
119            result = smatrixldltdet(ref a, ref pivots, n, isupper);
120            return result;
121        }
122
123
124        public static double determinantldlt(ref double[,] a,
125            ref int[] pivots,
126            int n,
127            bool isupper)
128        {
129            double result = 0;
130            int k = 0;
131
132            result = 1;
133            if( isupper )
134            {
135                k = 1;
136                while( k<=n )
137                {
138                    if( pivots[k]>0 )
139                    {
140                        result = result*a[k,k];
141                        k = k+1;
142                    }
143                    else
144                    {
145                        result = result*(a[k,k]*a[k+1,k+1]-a[k,k+1]*a[k,k+1]);
146                        k = k+2;
147                    }
148                }
149            }
150            else
151            {
152                k = n;
153                while( k>=1 )
154                {
155                    if( pivots[k]>0 )
156                    {
157                        result = result*a[k,k];
158                        k = k-1;
159                    }
160                    else
161                    {
162                        result = result*(a[k-1,k-1]*a[k,k]-a[k,k-1]*a[k,k-1]);
163                        k = k-2;
164                    }
165                }
166            }
167            return result;
168        }
169
170
171        public static double determinantsymmetric(double[,] a,
172            int n,
173            bool isupper)
174        {
175            double result = 0;
176            int[] pivots = new int[0];
177
178            a = (double[,])a.Clone();
179
180            ldlt.ldltdecomposition(ref a, n, isupper, ref pivots);
181            result = determinantldlt(ref a, ref pivots, n, isupper);
182            return result;
183        }
184    }
185}
Note: See TracBrowser for help on using the repository browser.