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
1314 slepc 1
/*
2
     SVD routines for accessing the problem matrix.
1376 slepc 3
 
4
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5
      SLEPc - Scalable Library for Eigenvalue Problem Computations
6
      Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
7
 
8
      This file is part of SLEPc. See the README file for conditions of use
9
      and additional information.
10
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1314 slepc 11
*/
1376 slepc 12
 
1521 slepc 13
#include "private/svdimpl.h"      /*I "slepcsvd.h" I*/
1314 slepc 14
 
15
#undef __FUNCT__
16
#define __FUNCT__ "SVDMatMult"
17
PetscErrorCode SVDMatMult(SVD svd,PetscTruth trans,Vec x,Vec y)
18
{
19
  PetscErrorCode ierr;
20
 
21
  PetscFunctionBegin;
22
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
23
  PetscValidHeaderSpecific(x,VEC_COOKIE,3);
24
  PetscValidHeaderSpecific(y,VEC_COOKIE,4);
25
  svd->matvecs++;
26
  if (trans) {
27
    if (svd->AT) {
28
      ierr = MatMult(svd->AT,x,y);CHKERRQ(ierr);
29
    } else {
30
      ierr = MatMultTranspose(svd->A,x,y);CHKERRQ(ierr);
31
    }
32
  } else {
33
    if (svd->A) {
34
      ierr = MatMult(svd->A,x,y);CHKERRQ(ierr);
35
    } else {
36
      ierr = MatMultTranspose(svd->AT,x,y);CHKERRQ(ierr);
37
    }
38
  }
39
  PetscFunctionReturn(0);
40
}
41
 
42
#undef __FUNCT__
43
#define __FUNCT__ "SVDMatGetVecs"
44
PetscErrorCode SVDMatGetVecs(SVD svd,Vec *x,Vec *y)
45
{
46
  PetscErrorCode ierr;
47
 
48
  PetscFunctionBegin;
49
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
50
  if (svd->A) {
51
    ierr = MatGetVecs(svd->A,x,y);CHKERRQ(ierr);
52
  } else {
53
    ierr = MatGetVecs(svd->AT,y,x);CHKERRQ(ierr);
54
  }
55
  PetscFunctionReturn(0);
56
}
57
 
58
#undef __FUNCT__
59
#define __FUNCT__ "SVDMatGetSize"
60
PetscErrorCode SVDMatGetSize(SVD svd,PetscInt *m,PetscInt *n)
61
{
62
  PetscErrorCode ierr;
63
 
64
  PetscFunctionBegin;
65
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
66
  if (svd->A) {
67
    ierr = MatGetSize(svd->A,m,n);CHKERRQ(ierr);
68
  } else {
69
    ierr = MatGetSize(svd->AT,n,m);CHKERRQ(ierr);
70
  }
71
  PetscFunctionReturn(0);
72
}
73
 
74
#undef __FUNCT__
75
#define __FUNCT__ "SVDMatGetLocalSize"
76
PetscErrorCode SVDMatGetLocalSize(SVD svd,PetscInt *m,PetscInt *n)
77
{
78
  PetscErrorCode ierr;
79
 
80
  PetscFunctionBegin;
81
  PetscValidHeaderSpecific(svd,SVD_COOKIE,1);
82
  if (svd->A) {
83
    ierr = MatGetLocalSize(svd->A,m,n);CHKERRQ(ierr);
84
  } else {
85
    ierr = MatGetLocalSize(svd->AT,n,m);CHKERRQ(ierr);
86
  }
87
  PetscFunctionReturn(0);
88
}