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