Subversion Repositories slepc-dev

Rev

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

Rev Author Line No. Line
1249 slepc 1
/*
2
     SVD routines for setting up the solver.
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
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1249 slepc 22
*/
1376 slepc 23
 
2284 jroman 24
#include <private/svdimpl.h>      /*I "slepcsvd.h" I*/
1249 slepc 25
 
26
#undef __FUNCT__  
27
#define __FUNCT__ "SVDSetOperator"
1333 slepc 28
/*@
1249 slepc 29
   SVDSetOperator - Set the matrix associated with the singular value problem.
30
 
31
   Collective on SVD and Mat
32
 
33
   Input Parameters:
34
+  svd - the singular value solver context
35
-  A  - the matrix associated with the singular value problem
36
 
37
   Level: beginner
38
 
1268 slepc 39
.seealso: SVDSolve(), SVDGetOperator()
1249 slepc 40
@*/
41
PetscErrorCode SVDSetOperator(SVD svd,Mat mat)
42
{
43
  PetscErrorCode ierr;
44
 
45
  PetscFunctionBegin;
2213 jroman 46
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
47
  PetscValidHeaderSpecific(mat,MAT_CLASSID,2);
1249 slepc 48
  PetscCheckSameComm(svd,1,mat,2);
2350 jroman 49
  if (svd->setupcalled) { ierr = SVDReset(svd);CHKERRQ(ierr); }
1249 slepc 50
  ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
2305 jroman 51
  ierr = MatDestroy(&svd->OP);CHKERRQ(ierr);
1314 slepc 52
  svd->OP = mat;
1249 slepc 53
  PetscFunctionReturn(0);
54
}
55
 
56
#undef __FUNCT__  
1268 slepc 57
#define __FUNCT__ "SVDGetOperator"
1333 slepc 58
/*@
59
   SVDGetOperator - Get the matrix associated with the singular value problem.
1255 slepc 60
 
1249 slepc 61
   Not collective, though parallel Mats are returned if the SVD is parallel
62
 
63
   Input Parameter:
64
.  svd - the singular value solver context
65
 
66
   Output Parameters:
1268 slepc 67
.  A    - the matrix associated with the singular value problem
1249 slepc 68
 
1255 slepc 69
   Level: advanced
1249 slepc 70
 
1268 slepc 71
.seealso: SVDSolve(), SVDSetOperator()
1249 slepc 72
@*/
1268 slepc 73
PetscErrorCode SVDGetOperator(SVD svd,Mat *A)
1249 slepc 74
{
75
  PetscFunctionBegin;
2213 jroman 76
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1268 slepc 77
  PetscValidPointer(A,2);
1314 slepc 78
  *A = svd->OP;
1249 slepc 79
  PetscFunctionReturn(0);
80
}
81
 
82
#undef __FUNCT__  
83
#define __FUNCT__ "SVDSetUp"
84
/*@
85
   SVDSetUp - Sets up all the internal data structures necessary for the
86
   execution of the singular value solver.
87
 
88
   Collective on SVD
89
 
90
   Input Parameter:
1405 slepc 91
.  svd   - singular value solver context
1249 slepc 92
 
93
   Level: advanced
94
 
95
   Notes:
96
   This function need not be called explicitly in most cases, since SVDSolve()
97
   calls it. It can be useful when one wants to measure the set-up time
98
   separately from the solve time.
99
 
100
.seealso: SVDCreate(), SVDSolve(), SVDDestroy()
101
@*/
102
PetscErrorCode SVDSetUp(SVD svd)
103
{
104
  PetscErrorCode ierr;
2216 jroman 105
  PetscBool      flg,lindep;
2410 jroman 106
  PetscInt       i,k,M,N;
1952 jroman 107
  PetscReal      norm;
1249 slepc 108
 
109
  PetscFunctionBegin;
2213 jroman 110
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1249 slepc 111
  if (svd->setupcalled) PetscFunctionReturn(0);
112
  ierr = PetscLogEventBegin(SVD_SetUp,svd,0,0,0);CHKERRQ(ierr);
2381 jroman 113
  if (!svd->ip) { ierr = SVDGetIP(svd,&svd->ip);CHKERRQ(ierr); }
1249 slepc 114
 
115
  /* Set default solver type */
1422 slepc 116
  if (!((PetscObject)svd)->type_name) {
1342 slepc 117
    ierr = SVDSetType(svd,SVDCROSS);CHKERRQ(ierr);
1249 slepc 118
  }
119
 
120
  /* check matrix */
2350 jroman 121
  if (!svd->OP) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_WRONGSTATE,"SVDSetOperator must be called first");
1256 slepc 122
 
1264 slepc 123
  /* determine how to build the transpose */
1283 slepc 124
  if (svd->transmode == PETSC_DECIDE) {
1314 slepc 125
    ierr = MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);CHKERRQ(ierr);    
1264 slepc 126
    if (flg) svd->transmode = SVD_TRANSPOSE_EXPLICIT;
1283 slepc 127
    else svd->transmode = SVD_TRANSPOSE_IMPLICIT;
1264 slepc 128
  }
129
 
1256 slepc 130
  /* build transpose matrix */
1314 slepc 131
  ierr = MatGetSize(svd->OP,&M,&N);CHKERRQ(ierr);
132
  ierr = PetscObjectReference((PetscObject)svd->OP);CHKERRQ(ierr);
1256 slepc 133
  switch (svd->transmode) {
134
    case SVD_TRANSPOSE_EXPLICIT:
1314 slepc 135
      ierr = MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);CHKERRQ(ierr);
2214 jroman 136
      if (!flg) SETERRQ(((PetscObject)svd)->comm,1,"Matrix has not defined the MatTranpose operation");
1314 slepc 137
      if (M>=N) {
138
        svd->A = svd->OP;
2331 jroman 139
        ierr = MatTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->AT);CHKERRQ(ierr);
1314 slepc 140
      } else {
2331 jroman 141
        ierr = MatTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->A);CHKERRQ(ierr);
1314 slepc 142
        svd->AT = svd->OP;
143
      }
1256 slepc 144
      break;
1283 slepc 145
    case SVD_TRANSPOSE_IMPLICIT:
1314 slepc 146
      ierr = MatHasOperation(svd->OP,MATOP_MULT_TRANSPOSE,&flg);CHKERRQ(ierr);
2214 jroman 147
      if (!flg) SETERRQ(((PetscObject)svd)->comm,1,"Matrix has not defined the MatMultTranpose operation");
1314 slepc 148
      if (M>=N) {
149
        svd->A = svd->OP;
150
        svd->AT = PETSC_NULL;    
151
      } else {
152
        svd->A = PETSC_NULL;
153
        svd->AT = svd->OP;
154
      }
1256 slepc 155
      break;
156
    default:
2214 jroman 157
      SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid transpose mode");
1256 slepc 158
  }
1249 slepc 159
 
2410 jroman 160
  if (svd->A) {
161
    ierr = SlepcMatGetVecsTemplate(svd->A,&svd->tr,&svd->tl);CHKERRQ(ierr);
162
  } else {
163
    ierr = SlepcMatGetVecsTemplate(svd->AT,&svd->tl,&svd->tr);CHKERRQ(ierr);
164
  }
165
 
1249 slepc 166
  /* call specific solver setup */
167
  ierr = (*svd->ops->setup)(svd);CHKERRQ(ierr);
168
 
1272 slepc 169
  if (svd->ncv > M || svd->ncv > N)
2214 jroman 170
    SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"ncv bigger than matrix dimensions");
1272 slepc 171
  if (svd->nsv > svd->ncv)
2214 jroman 172
    SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nsv bigger than ncv");
1272 slepc 173
 
174
  if (svd->ncv != svd->n) {
175
    /* free memory for previous solution  */
176
    if (svd->n) {
177
      ierr = PetscFree(svd->sigma);CHKERRQ(ierr);
1603 slepc 178
      ierr = PetscFree(svd->perm);CHKERRQ(ierr);
1288 slepc 179
      ierr = PetscFree(svd->errest);CHKERRQ(ierr);
2410 jroman 180
      ierr = VecDestroyVecs(svd->n,&svd->V);CHKERRQ(ierr);
1272 slepc 181
    }
182
    /* allocate memory for next solution */
183
    ierr = PetscMalloc(svd->ncv*sizeof(PetscReal),&svd->sigma);CHKERRQ(ierr);
1889 jroman 184
    ierr = PetscMalloc(svd->ncv*sizeof(PetscInt),&svd->perm);CHKERRQ(ierr);
1288 slepc 185
    ierr = PetscMalloc(svd->ncv*sizeof(PetscReal),&svd->errest);CHKERRQ(ierr);
2410 jroman 186
    ierr = VecDuplicateVecs(svd->tr,svd->ncv,&svd->V);CHKERRQ(ierr);
1272 slepc 187
    svd->n = svd->ncv;
1263 slepc 188
  }
189
 
1952 jroman 190
  /* process initial vectors */
191
  if (svd->nini<0) {
192
    svd->nini = -svd->nini;
2219 jroman 193
    if (svd->nini>svd->ncv) SETERRQ(((PetscObject)svd)->comm,1,"The number of initial vectors is larger than ncv");
1952 jroman 194
    k = 0;
195
    for (i=0;i<svd->nini;i++) {
196
      ierr = VecCopy(svd->IS[i],svd->V[k]);CHKERRQ(ierr);
2305 jroman 197
      ierr = VecDestroy(&svd->IS[i]);CHKERRQ(ierr);
1952 jroman 198
      ierr = IPOrthogonalize(svd->ip,0,PETSC_NULL,k,PETSC_NULL,svd->V,svd->V[k],PETSC_NULL,&norm,&lindep);CHKERRQ(ierr);
199
      if (norm==0.0 || lindep) PetscInfo(svd,"Linearly dependent initial vector found, removing...\n");
200
      else {
201
        ierr = VecScale(svd->V[k],1.0/norm);CHKERRQ(ierr);
202
        k++;
203
      }
204
    }
205
    svd->nini = k;
206
    ierr = PetscFree(svd->IS);CHKERRQ(ierr);
207
  }
208
 
1249 slepc 209
  ierr = PetscLogEventEnd(SVD_SetUp,svd,0,0,0);CHKERRQ(ierr);
210
  svd->setupcalled = 1;
211
  PetscFunctionReturn(0);
212
}
1952 jroman 213
 
214
#undef __FUNCT__  
215
#define __FUNCT__ "SVDSetInitialSpace"
216
/*@
217
   SVDSetInitialSpace - Specify a basis of vectors that constitute the initial
218
   space, that is, the subspace from which the solver starts to iterate.
219
 
220
   Collective on SVD and Vec
221
 
222
   Input Parameter:
223
+  svd   - the singular value solver context
224
.  n     - number of vectors
225
-  is    - set of basis vectors of the initial space
226
 
227
   Notes:
228
   Some solvers start to iterate on a single vector (initial vector). In that case,
229
   the other vectors are ignored.
230
 
231
   These vectors do not persist from one SVDSolve() call to the other, so the
232
   initial space should be set every time.
233
 
234
   The vectors do not need to be mutually orthonormal, since they are explicitly
235
   orthonormalized internally.
236
 
237
   Common usage of this function is when the user can provide a rough approximation
238
   of the wanted singular space. Then, convergence may be faster.
239
 
240
   Level: intermediate
241
@*/
242
PetscErrorCode SVDSetInitialSpace(SVD svd,PetscInt n,Vec *is)
243
{
244
  PetscErrorCode ierr;
245
  PetscInt       i;
246
 
247
  PetscFunctionBegin;
2213 jroman 248
  PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
2326 jroman 249
  PetscValidLogicalCollectiveInt(svd,n,2);
2214 jroman 250
  if (n<0) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
1952 jroman 251
 
252
  /* free previous non-processed vectors */
253
  if (svd->nini<0) {
254
    for (i=0;i<-svd->nini;i++) {
2305 jroman 255
      ierr = VecDestroy(&svd->IS[i]);CHKERRQ(ierr);
1952 jroman 256
    }
257
    ierr = PetscFree(svd->IS);CHKERRQ(ierr);
258
  }
259
 
260
  /* get references of passed vectors */
261
  ierr = PetscMalloc(n*sizeof(Vec),&svd->IS);CHKERRQ(ierr);
262
  for (i=0;i<n;i++) {
263
    ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);
264
    svd->IS[i] = is[i];
265
  }
266
 
267
  svd->nini = -n;
268
  svd->setupcalled = 0;
269
  PetscFunctionReturn(0);
270
}
271