Welcome to the NetCologne GmbH open source mirroring service!

This machine mirrors various open-source projects. 20 Gbit/s uplink.

If there are any issues or you want another project mirrored, please contact mirror-service -=AT=- netcologne DOT de !

GetFEM: src/gmm/gmm_dense_lu.h Source File
GetFEM  5.4.2
gmm_dense_lu.h
Go to the documentation of this file.
1 /* -*- c++ -*- (enables emacs c++ mode) */
2 /*===========================================================================
3 
4  Copyright (C) 2003-2020 Yves Renard
5 
6  This file is a part of GetFEM
7 
8  GetFEM is free software; you can redistribute it and/or modify it
9  under the terms of the GNU Lesser General Public License as published
10  by the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version along with the GCC Runtime Library
12  Exception either version 3.1 or (at your option) any later version.
13  This program is distributed in the hope that it will be useful, but
14  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
15  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
16  License and GCC Runtime Library Exception for more details.
17  You should have received a copy of the GNU Lesser General Public License
18  along with this program; if not, write to the Free Software Foundation,
19  Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.
20 
21  As a special exception, you may use this file as it is a part of a free
22  software library without restriction. Specifically, if other files
23  instantiate templates or use macros or inline functions from this file,
24  or you compile this file and link it with other files to produce an
25  executable, this file does not by itself cause the resulting executable
26  to be covered by the GNU Lesser General Public License. This exception
27  does not however invalidate any other reasons why the executable file
28  might be covered by the GNU Lesser General Public License.
29 
30 ===========================================================================*/
31 
32 // This file is a modified version of lu.h from MTL.
33 // See http://osl.iu.edu/research/mtl/
34 // Following the corresponding Copyright notice.
35 //===========================================================================
36 //
37 // Copyright (c) 1998-2020, University of Notre Dame. All rights reserved.
38 // Redistribution and use in source and binary forms, with or without
39 // modification, are permitted provided that the following conditions are met:
40 //
41 // * Redistributions of source code must retain the above copyright
42 // notice, this list of conditions and the following disclaimer.
43 // * Redistributions in binary form must reproduce the above copyright
44 // notice, this list of conditions and the following disclaimer in the
45 // documentation and/or other materials provided with the distribution.
46 // * Neither the name of the University of Notre Dame nor the
47 // names of its contributors may be used to endorse or promote products
48 // derived from this software without specific prior written permission.
49 //
50 // THIS SOFTWARE IS PROVIDED BY THE TRUSTEES OF INDIANA UNIVERSITY AND
51 // CONTRIBUTORS ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING,
52 // BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
53 // FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE TRUSTEES
54 // OF INDIANA UNIVERSITY AND CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
55 // INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT
56 // NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
57 // DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
58 // THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
59 // (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF
60 // THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
61 //
62 //===========================================================================
63 
64 /**@file gmm_dense_lu.h
65  @author Andrew Lumsdaine, Jeremy G. Siek, Lie-Quan Lee, Y. Renard
66  @date June 5, 2003.
67  @brief LU factorizations and determinant computation for dense matrices.
68 */
69 #ifndef GMM_DENSE_LU_H
70 #define GMM_DENSE_LU_H
71 
72 #include "gmm_dense_Householder.h"
73 
74 namespace gmm {
75 
76  /* ********************************************************************** */
77  /* IPVT structure. */
78  /* ********************************************************************** */
79  // For compatibility with lapack version with 64 or 32 bit integer.
80  // Should be replaced by std::vector<size_type> if 32 bit integer version
81  // of lapack is not used anymore (and lapack_ipvt_int set to size_type)
82 
83  // Do not use iterators of this interface container
84  class lapack_ipvt : public std::vector<size_type> {
85  bool is_int64;
86  size_type &operator[](size_type i)
87  { return std::vector<size_type>::operator[](i); }
88  size_type operator[] (size_type i) const
89  { return std::vector<size_type>::operator[](i); }
90  void begin(void) const {}
91  void begin(void) {}
92  void end(void) const {}
93  void end(void) {}
94 
95  public:
96  void set_to_int32() { is_int64 = false; }
97  const size_type *pfirst() const
98  { return &(*(std::vector<size_type>::begin())); }
99  size_type *pfirst() { return &(*(std::vector<size_type>::begin())); }
100 
101  lapack_ipvt(size_type n) : std::vector<size_type>(n), is_int64(true) {}
102 
103  size_type get(size_type i) const {
104  const size_type *p = pfirst();
105  return is_int64 ? p[i] : size_type(((const int *)(p))[i]);
106  }
107  void set(size_type i, size_type val) {
108  size_type *p = pfirst();
109  if (is_int64) p[i] = val; else ((int *)(p))[i] = int(val);
110  }
111  };
112 }
113 
114 #include "gmm_opt.h"
115 
116 namespace gmm {
117 
118  /** LU Factorization of a general (dense) matrix (real or complex).
119 
120  This is the outer product (a level-2 operation) form of the LU
121  Factorization with pivoting algorithm . This is equivalent to
122  LAPACK's dgetf2. Also see "Matrix Computations" 3rd Ed. by Golub
123  and Van Loan section 3.2.5 and especially page 115.
124 
125  The pivot indices in ipvt are indexed starting from 1
126  so that this is compatible with LAPACK (Fortran).
127  */
128  template <typename DenseMatrix>
129  size_type lu_factor(DenseMatrix& A, lapack_ipvt& ipvt) {
130  typedef typename linalg_traits<DenseMatrix>::value_type T;
131  typedef typename number_traits<T>::magnitude_type R;
132  size_type info(0), i, j, jp, M(mat_nrows(A)), N(mat_ncols(A));
133  size_type NN = std::min(M, N);
134  std::vector<T> c(M), r(N);
135 
136  GMM_ASSERT2(ipvt.size()+1 >= NN, "IPVT too small");
137  for (i = 0; i+1 < NN; ++i) ipvt.set(i, i);
138 
139  if (M || N) {
140  for (j = 0; j+1 < NN; ++j) {
141  R max = gmm::abs(A(j,j)); jp = j;
142  for (i = j+1; i < M; ++i) /* find pivot. */
143  if (gmm::abs(A(i,j)) > max) { jp = i; max = gmm::abs(A(i,j)); }
144  ipvt.set(j, jp + 1);
145 
146  if (max == R(0)) { info = j + 1; break; }
147  if (jp != j) for (i = 0; i < N; ++i) std::swap(A(jp, i), A(j, i));
148 
149  for (i = j+1; i < M; ++i) { A(i, j) /= A(j,j); c[i-j-1] = -A(i, j); }
150  for (i = j+1; i < N; ++i) r[i-j-1] = A(j, i); // avoid the copy ?
151  rank_one_update(sub_matrix(A, sub_interval(j+1, M-j-1),
152  sub_interval(j+1, N-j-1)), c, conjugated(r));
153  }
154  ipvt.set(NN-1, NN);
155  }
156  return info;
157  }
158 
159  /** LU Solve : Solve equation Ax=b, given an LU factored matrix.*/
160  // Thanks to Valient Gough for this routine!
161  template <typename DenseMatrix, typename VectorB, typename VectorX,
162  typename Pvector>
163  void lu_solve(const DenseMatrix &LU, const Pvector& pvector,
164  VectorX &x, const VectorB &b) {
165  typedef typename linalg_traits<DenseMatrix>::value_type T;
166  copy(b, x);
167  for(size_type i = 0; i < pvector.size(); ++i) {
168  size_type perm = pvector.get(i)-1; // permutations stored in 1's offset
169  if(i != perm) { T aux = x[i]; x[i] = x[perm]; x[perm] = aux; }
170  }
171  /* solve Ax = b -> LUx = b -> Ux = L^-1 b. */
172  lower_tri_solve(LU, x, true);
173  upper_tri_solve(LU, x, false);
174  }
175 
176  template <typename DenseMatrix, typename VectorB, typename VectorX>
177  void lu_solve(const DenseMatrix &A, VectorX &x, const VectorB &b) {
178  typedef typename linalg_traits<DenseMatrix>::value_type T;
179  dense_matrix<T> B(mat_nrows(A), mat_ncols(A));
180  lapack_ipvt ipvt(mat_nrows(A));
181  gmm::copy(A, B);
182  size_type info = lu_factor(B, ipvt);
183  GMM_ASSERT1(!info, "Singular system, pivot = " << info);
184  lu_solve(B, ipvt, x, b);
185  }
186 
187  template <typename DenseMatrix, typename VectorB, typename VectorX,
188  typename Pvector>
189  void lu_solve_transposed(const DenseMatrix &LU, const Pvector& pvector,
190  VectorX &x, const VectorB &b) {
191  typedef typename linalg_traits<DenseMatrix>::value_type T;
192  copy(b, x);
193  lower_tri_solve(transposed(LU), x, false);
194  upper_tri_solve(transposed(LU), x, true);
195  for(size_type i = pvector.size(); i > 0; --i) {
196  size_type perm = pvector.get(i-1)-1; // permutations stored in 1's offset
197  if(i-1 != perm) { T aux = x[i-1]; x[i-1] = x[perm]; x[perm] = aux; }
198  }
199  }
200 
201 
202  ///@cond DOXY_SHOW_ALL_FUNCTIONS
203  template <typename DenseMatrixLU, typename DenseMatrix, typename Pvector>
204  void lu_inverse(const DenseMatrixLU& LU, const Pvector& pvector,
205  DenseMatrix& AInv, col_major) {
206  typedef typename linalg_traits<DenseMatrixLU>::value_type T;
207  std::vector<T> tmp(pvector.size(), T(0));
208  std::vector<T> result(pvector.size());
209  for(size_type i = 0; i < pvector.size(); ++i) {
210  tmp[i] = T(1);
211  lu_solve(LU, pvector, result, tmp);
212  copy(result, mat_col(AInv, i));
213  tmp[i] = T(0);
214  }
215  }
216 
217  template <typename DenseMatrixLU, typename DenseMatrix, typename Pvector>
218  void lu_inverse(const DenseMatrixLU& LU, const Pvector& pvector,
219  DenseMatrix& AInv, row_major) {
220  typedef typename linalg_traits<DenseMatrixLU>::value_type T;
221  std::vector<T> tmp(pvector.size(), T(0));
222  std::vector<T> result(pvector.size());
223  for(size_type i = 0; i < pvector.size(); ++i) {
224  tmp[i] = T(1); // to be optimized !!
225  // on peut sur le premier tri solve reduire le systeme
226  // et peut etre faire un solve sur une serie de vecteurs au lieu
227  // de vecteur a vecteur (accumulation directe de l'inverse dans la
228  // matrice au fur et a mesure du calcul ... -> evite la copie finale
229  lu_solve_transposed(LU, pvector, result, tmp);
230  copy(result, mat_row(AInv, i));
231  tmp[i] = T(0);
232  }
233  }
234  ///@endcond
235 
236  /** Given an LU factored matrix, build the inverse of the matrix. */
237  template <typename DenseMatrixLU, typename DenseMatrix, typename Pvector>
238  void lu_inverse(const DenseMatrixLU& LU, const Pvector& pvector,
239  const DenseMatrix& AInv_) {
240  DenseMatrix& AInv = const_cast<DenseMatrix&>(AInv_);
241  lu_inverse(LU, pvector, AInv, typename principal_orientation_type<typename
242  linalg_traits<DenseMatrix>::sub_orientation>::potype());
243  }
244 
245  /** Given a dense matrix, build the inverse of the matrix, and
246  return the determinant */
247  template <typename DenseMatrix>
248  typename linalg_traits<DenseMatrix>::value_type
249  lu_inverse(const DenseMatrix& A_, bool doassert = true) {
250  typedef typename linalg_traits<DenseMatrix>::value_type T;
251  DenseMatrix& A = const_cast<DenseMatrix&>(A_);
252  dense_matrix<T> B(mat_nrows(A), mat_ncols(A));
253  lapack_ipvt ipvt(mat_nrows(A));
254  gmm::copy(A, B);
255  size_type info = lu_factor(B, ipvt);
256  if (doassert) GMM_ASSERT1(!info, "Non invertible matrix, pivot = "<<info);
257  if (!info) lu_inverse(B, ipvt, A);
258  return lu_det(B, ipvt);
259  }
260 
261  /** Compute the matrix determinant (via a LU factorization) */
262  template <typename DenseMatrixLU, typename Pvector>
263  typename linalg_traits<DenseMatrixLU>::value_type
264  lu_det(const DenseMatrixLU& LU, const Pvector &pvector) {
265  typedef typename linalg_traits<DenseMatrixLU>::value_type T;
266  T det(1);
267  for (size_type j = 0; j < std::min(mat_nrows(LU), mat_ncols(LU)); ++j)
268  det *= LU(j,j);
269  for(size_type i = 0; i < pvector.size(); ++i)
270  if (i != size_type(pvector.get(i)-1)) { det = -det; }
271  return det;
272  }
273 
274  template <typename DenseMatrix>
275  typename linalg_traits<DenseMatrix>::value_type
276  lu_det(const DenseMatrix& A) {
277  typedef typename linalg_traits<DenseMatrix>::value_type T;
278  dense_matrix<T> B(mat_nrows(A), mat_ncols(A));
279  lapack_ipvt ipvt(mat_nrows(A));
280  gmm::copy(A, B);
281  lu_factor(B, ipvt);
282  return lu_det(B, ipvt);
283  }
284 
285 }
286 
287 #endif
288 
bgeot::size_type
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:49
gmm::conjugated
conjugated_return< L >::return_type conjugated(const L &v)
return a conjugated view of the input matrix or vector.
Definition: gmm_conjugated.h:294
gmm_dense_Householder.h
Householder for dense matrices.
gmm_opt.h
Optimization for some small cases (inversion of 2x2 matrices etc.)
gmm::lu_solve
void lu_solve(const DenseMatrix &LU, const Pvector &pvector, VectorX &x, const VectorB &b)
LU Solve : Solve equation Ax=b, given an LU factored matrix.
Definition: gmm_dense_lu.h:163
gmm::copy
void copy(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:977