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