Subversion Repositories slepc-dev

Rev

Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1298 slepc 1
/*                      
2
     This file contains routines for handling small-size dense problems.
3
     All routines are simply wrappers to LAPACK routines.
1376 slepc 4
 
5
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 6
   SLEPc - Scalable Library for Eigenvalue Problem Computations
7
   Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
1376 slepc 8
 
1672 slepc 9
   This file is part of SLEPc.
10
 
11
   SLEPc is free software: you can redistribute it and/or modify it under  the
12
   terms of version 3 of the GNU Lesser General Public License as published by
13
   the Free Software Foundation.
14
 
15
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
16
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
17
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
18
   more details.
19
 
20
   You  should have received a copy of the GNU Lesser General  Public  License
21
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
1376 slepc 22
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1298 slepc 23
*/
24
 
1521 slepc 25
#include "private/svdimpl.h"        /*I "slepcsvd.h" I*/
1298 slepc 26
#include "slepcblaslapack.h"
27
 
1333 slepc 28
#undef __FUNCT__  
29
#define __FUNCT__ "SVDDense"
30
/*@
31
   SVDDense - Solves a dense singular value problem.
32
 
33
   Not Collective
34
 
35
   Input Parameters:
36
+  M  - dimension of the problem (rows)
37
.  N  - dimension of the problem (colums)
38
-  A  - pointer to the array containing the matrix values
39
 
40
   Output Parameters:
41
+  sigma  - pointer to the array to store the computed singular values
42
.  U  - pointer to the array to store left singular vectors
43
-  VT  - pointer to the array to store right singular vectors
44
 
45
   Matrix A is overwritten.
46
 
47
   This routine uses LAPACK routines xGESDD.
48
 
49
   Level: developer
50
 
51
@*/
1509 slepc 52
PetscErrorCode SVDDense(PetscInt M_,PetscInt N_,PetscScalar* A,PetscReal* sigma,PetscScalar* U,PetscScalar* VT)
1298 slepc 53
{
1336 slepc 54
#if defined(SLEPC_MISSING_LAPACK_GESDD)
55
  PetscFunctionBegin;
56
  SETERRQ(PETSC_ERR_SUP,"GESDD - Lapack routine is unavailable.");
57
#else
1298 slepc 58
  PetscErrorCode ierr;
59
  PetscScalar    qwork,*work;
1643 slepc 60
  PetscBLASInt   n,info,lwork,*iwork,M,N;
1298 slepc 61
#if defined(PETSC_USE_COMPLEX)
62
  PetscReal       *rwork;
63
#endif
64
 
65
  PetscFunctionBegin;
66
  /* workspace query & allocation */
1339 slepc 67
  ierr = PetscLogEventBegin(SVD_Dense,0,0,0,0);CHKERRQ(ierr);
1643 slepc 68
  M = PetscBLASIntCast(M_);
69
  N = PetscBLASIntCast(N_);
1298 slepc 70
  n = PetscMin(M,N);
1504 slepc 71
  ierr = PetscMalloc(sizeof(PetscInt)*8*n,&iwork);CHKERRQ(ierr);
1298 slepc 72
  lwork = -1;
73
#if defined(PETSC_USE_COMPLEX)
74
  ierr = PetscMalloc(sizeof(PetscReal)*(5*n*n+7*n),&rwork);CHKERRQ(ierr);
1536 slepc 75
  LAPACKgesdd_("O",&M,&N,A,&M,sigma,U,&M,VT,&N,&qwork,&lwork,rwork,iwork,&info);
1298 slepc 76
#else
1536 slepc 77
  LAPACKgesdd_("O",&M,&N,A,&M,sigma,U,&M,VT,&N,&qwork,&lwork,iwork,&info);
1298 slepc 78
#endif
79
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGESDD %d",info);
1643 slepc 80
  lwork = (PetscBLASInt)PetscRealPart(qwork);
1298 slepc 81
  ierr = PetscMalloc(sizeof(PetscScalar)*lwork,&work);CHKERRQ(ierr);
82
 
83
  /* computation */  
84
#if defined(PETSC_USE_COMPLEX)
1536 slepc 85
  LAPACKgesdd_("O",&M,&N,A,&M,sigma,U,&M,VT,&N,work,&lwork,rwork,iwork,&info);
1298 slepc 86
  ierr = PetscFree(rwork);CHKERRQ(ierr);
87
#else
1536 slepc 88
  LAPACKgesdd_("O",&M,&N,A,&M,sigma,U,&M,VT,&N,work,&lwork,iwork,&info);
1298 slepc 89
#endif
90
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGESDD %d",info);
91
  ierr = PetscFree(iwork);CHKERRQ(ierr);
92
  ierr = PetscFree(work);CHKERRQ(ierr);
1339 slepc 93
  ierr = PetscLogEventEnd(SVD_Dense,0,0,0,0);CHKERRQ(ierr);
1298 slepc 94
  PetscFunctionReturn(0);
1336 slepc 95
#endif
1326 slepc 96
}