Subversion Repositories slepc-dev

Rev

Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1314 slepc 1
/*
2
     SVD routines for accessing the problem matrix.
1376 slepc 3
 
4
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 5
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2116 eromero 6
   Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
1376 slepc 7
 
1672 slepc 8
   This file is part of SLEPc.
9
 
10
   SLEPc is free software: you can redistribute it and/or modify it under  the
11
   terms of version 3 of the GNU Lesser General Public License as published by
12
   the Free Software Foundation.
13
 
14
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
15
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
16
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
17
   more details.
18
 
19
   You  should have received a copy of the GNU Lesser General  Public  License
20
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
1376 slepc 21
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1314 slepc 22
*/
1376 slepc 23
 
1521 slepc 24
#include "private/svdimpl.h"      /*I "slepcsvd.h" I*/
1314 slepc 25
 
26
#undef __FUNCT__
27
#define __FUNCT__ "SVDMatMult"
28
PetscErrorCode SVDMatMult(SVD svd,PetscTruth trans,Vec x,Vec y)
29
{
30
  PetscErrorCode ierr;
31
 
32
  PetscFunctionBegin;
33
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
34
  PetscValidHeaderSpecific(x,VEC_COOKIE,3);
35
  PetscValidHeaderSpecific(y,VEC_COOKIE,4);
36
  svd->matvecs++;
37
  if (trans) {
38
    if (svd->AT) {
39
      ierr = MatMult(svd->AT,x,y);CHKERRQ(ierr);
40
    } else {
41
      ierr = MatMultTranspose(svd->A,x,y);CHKERRQ(ierr);
42
    }
43
  } else {
44
    if (svd->A) {
45
      ierr = MatMult(svd->A,x,y);CHKERRQ(ierr);
46
    } else {
47
      ierr = MatMultTranspose(svd->AT,x,y);CHKERRQ(ierr);
48
    }
49
  }
50
  PetscFunctionReturn(0);
51
}
52
 
53
#undef __FUNCT__
54
#define __FUNCT__ "SVDMatGetVecs"
55
PetscErrorCode SVDMatGetVecs(SVD svd,Vec *x,Vec *y)
56
{
57
  PetscErrorCode ierr;
58
 
59
  PetscFunctionBegin;
60
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
61
  if (svd->A) {
62
    ierr = MatGetVecs(svd->A,x,y);CHKERRQ(ierr);
63
  } else {
64
    ierr = MatGetVecs(svd->AT,y,x);CHKERRQ(ierr);
65
  }
66
  PetscFunctionReturn(0);
67
}
68
 
69
#undef __FUNCT__
70
#define __FUNCT__ "SVDMatGetSize"
71
PetscErrorCode SVDMatGetSize(SVD svd,PetscInt *m,PetscInt *n)
72
{
73
  PetscErrorCode ierr;
74
 
75
  PetscFunctionBegin;
76
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
77
  if (svd->A) {
78
    ierr = MatGetSize(svd->A,m,n);CHKERRQ(ierr);
79
  } else {
80
    ierr = MatGetSize(svd->AT,n,m);CHKERRQ(ierr);
81
  }
82
  PetscFunctionReturn(0);
83
}
84
 
85
#undef __FUNCT__
86
#define __FUNCT__ "SVDMatGetLocalSize"
87
PetscErrorCode SVDMatGetLocalSize(SVD svd,PetscInt *m,PetscInt *n)
88
{
89
  PetscErrorCode ierr;
90
 
91
  PetscFunctionBegin;
92
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
93
  if (svd->A) {
94
    ierr = MatGetLocalSize(svd->A,m,n);CHKERRQ(ierr);
95
  } else {
96
    ierr = MatGetLocalSize(svd->AT,n,m);CHKERRQ(ierr);
97
  }
98
  PetscFunctionReturn(0);
99
}