Actual source code: svdsolve.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:       SVD routines related to the solution process.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 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.

 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.

 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/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <slepc/private/svdimpl.h>   /*I "slepcsvd.h" I*/

 28: PetscErrorCode SVDComputeVectors(SVD svd)
 29: {
 31:   Vec            tl,uj,vj;
 32:   PetscInt       j,oldsize;
 33:   PetscReal      norm;

 36:   SVDCheckSolved(svd,1);
 37:   switch (svd->state) {
 38:   case SVD_STATE_SOLVED:
 39:     /* generate left singular vectors on U */
 40:     if (!svd->U) { SVDGetBV(svd,NULL,&svd->U); }
 41:     BVGetSizes(svd->U,NULL,NULL,&oldsize);
 42:     if (!oldsize) {
 43:       SVDMatCreateVecs(svd,NULL,&tl);
 44:       BVSetSizesFromVec(svd->U,tl,svd->ncv);
 45:       VecDestroy(&tl);
 46:     }
 47:     for (j=0;j<svd->nconv;j++) {
 48:       BVGetColumn(svd->V,j,&vj);
 49:       BVGetColumn(svd->U,j,&uj);
 50:       SVDMatMult(svd,PETSC_FALSE,vj,uj);
 51:       BVRestoreColumn(svd->V,j,&vj);
 52:       BVRestoreColumn(svd->U,j,&uj);
 53:       BVOrthogonalizeColumn(svd->U,j,NULL,&norm,NULL);
 54:       BVScaleColumn(svd->U,j,1.0/norm);
 55:     }
 56:     break;
 57:   default:
 58:     break;
 59:   }
 60:   svd->state = SVD_STATE_VECTORS;
 61:   return(0);
 62: }

 66: /*@
 67:    SVDSolve - Solves the singular value problem.

 69:    Collective on SVD

 71:    Input Parameter:
 72: .  svd - singular value solver context obtained from SVDCreate()

 74:    Options Database Keys:
 75: +  -svd_view - print information about the solver used
 76: .  -svd_view_mat binary - save the matrix to the default binary viewer
 77: .  -svd_view_vectors binary - save the computed singular vectors to the default binary viewer
 78: .  -svd_view_values - print computed singular values
 79: .  -svd_converged_reason - print reason for convergence, and number of iterations
 80: .  -svd_error_absolute - print absolute errors of each singular triplet
 81: -  -svd_error_relative - print relative errors of each singular triplet

 83:    Level: beginner

 85: .seealso: SVDCreate(), SVDSetUp(), SVDDestroy()
 86: @*/
 87: PetscErrorCode SVDSolve(SVD svd)
 88: {
 90:   PetscInt       i,*workperm;

 94:   PetscLogEventBegin(SVD_Solve,svd,0,0,0);

 96:   /* call setup */
 97:   SVDSetUp(svd);
 98:   svd->its = 0;
 99:   svd->nconv = 0;
100:   for (i=0;i<svd->ncv;i++) {
101:     svd->sigma[i]  = 0.0;
102:     svd->errest[i] = 0.0;
103:     svd->perm[i]   = i;
104:   }
105:   SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,svd->ncv);
106:   SVDViewFromOptions(svd,NULL,"-svd_view_pre");

108:   (*svd->ops->solve)(svd);
109:   svd->state = (svd->leftbasis)? SVD_STATE_VECTORS: SVD_STATE_SOLVED;

111:   /* sort singular triplets */
112:   if (svd->which == SVD_SMALLEST) {
113:     PetscSortRealWithPermutation(svd->nconv,svd->sigma,svd->perm);
114:   } else {
115:     PetscMalloc(sizeof(PetscInt)*svd->nconv,&workperm);
116:     for (i=0;i<svd->nconv;i++) workperm[i] = i;
117:     PetscSortRealWithPermutation(svd->nconv,svd->sigma,workperm);
118:     for (i=0;i<svd->nconv;i++) svd->perm[i] = workperm[svd->nconv-i-1];
119:     PetscFree(workperm);
120:   }
121:   PetscLogEventEnd(SVD_Solve,svd,0,0,0);

123:   /* various viewers */
124:   SVDViewFromOptions(svd,NULL,"-svd_view");
125:   SVDReasonViewFromOptions(svd);
126:   SVDErrorViewFromOptions(svd);
127:   SVDValuesViewFromOptions(svd);
128:   SVDVectorsViewFromOptions(svd);
129:   MatViewFromOptions(svd->OP,(PetscObject)svd,"-svd_view_mat");

131:   /* Remove the initial subspaces */
132:   svd->nini = 0;
133:   svd->ninil = 0;
134:   return(0);
135: }

139: /*@
140:    SVDGetIterationNumber - Gets the current iteration number. If the
141:    call to SVDSolve() is complete, then it returns the number of iterations
142:    carried out by the solution method.

144:    Not Collective

146:    Input Parameter:
147: .  svd - the singular value solver context

149:    Output Parameter:
150: .  its - number of iterations

152:    Level: intermediate

154:    Note:
155:    During the i-th iteration this call returns i-1. If SVDSolve() is
156:    complete, then parameter "its" contains either the iteration number at
157:    which convergence was successfully reached, or failure was detected.
158:    Call SVDGetConvergedReason() to determine if the solver converged or
159:    failed and why.

161: .seealso: SVDGetConvergedReason(), SVDSetTolerances()
162: @*/
163: PetscErrorCode SVDGetIterationNumber(SVD svd,PetscInt *its)
164: {
168:   *its = svd->its;
169:   return(0);
170: }

174: /*@
175:    SVDGetConvergedReason - Gets the reason why the SVDSolve() iteration was
176:    stopped.

178:    Not Collective

180:    Input Parameter:
181: .  svd - the singular value solver context

183:    Output Parameter:
184: .  reason - negative value indicates diverged, positive value converged
185:    (see SVDConvergedReason)

187:    Possible values for reason:
188: +  SVD_CONVERGED_TOL - converged up to tolerance
189: .  SVD_DIVERGED_ITS - required more than its to reach convergence
190: -  SVD_DIVERGED_BREAKDOWN - generic breakdown in method

192:    Level: intermediate

194:    Notes: Can only be called after the call to SVDSolve() is complete.

196: .seealso: SVDSetTolerances(), SVDSolve(), SVDConvergedReason
197: @*/
198: PetscErrorCode SVDGetConvergedReason(SVD svd,SVDConvergedReason *reason)
199: {
203:   SVDCheckSolved(svd,1);
204:   *reason = svd->reason;
205:   return(0);
206: }

210: /*@
211:    SVDGetConverged - Gets the number of converged singular values.

213:    Not Collective

215:    Input Parameter:
216: .  svd - the singular value solver context

218:    Output Parameter:
219: .  nconv - number of converged singular values

221:    Note:
222:    This function should be called after SVDSolve() has finished.

224:    Level: beginner

226: @*/
227: PetscErrorCode SVDGetConverged(SVD svd,PetscInt *nconv)
228: {
232:   SVDCheckSolved(svd,1);
233:   *nconv = svd->nconv;
234:   return(0);
235: }

239: /*@
240:    SVDGetSingularTriplet - Gets the i-th triplet of the singular value decomposition
241:    as computed by SVDSolve(). The solution consists in the singular value and its left
242:    and right singular vectors.

244:    Not Collective, but vectors are shared by all processors that share the SVD

246:    Input Parameters:
247: +  svd - singular value solver context
248: -  i   - index of the solution

250:    Output Parameters:
251: +  sigma - singular value
252: .  u     - left singular vector
253: -  v     - right singular vector

255:    Note:
256:    Both U or V can be NULL if singular vectors are not required.
257:    Otherwise, the caller must provide valid Vec objects, i.e.,
258:    they must be created by the calling program with e.g. MatCreateVecs().

260:    The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
261:    Singular triplets are indexed according to the ordering criterion established
262:    with SVDSetWhichSingularTriplets().

264:    Level: beginner

266: .seealso: SVDSolve(), SVDGetConverged(), SVDSetWhichSingularTriplets()
267: @*/
268: PetscErrorCode SVDGetSingularTriplet(SVD svd,PetscInt i,PetscReal *sigma,Vec u,Vec v)
269: {
271:   PetscInt       M,N;
272:   Vec            w;

277:   SVDCheckSolved(svd,1);
280:   if (i<0 || i>=svd->nconv) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
281:   *sigma = svd->sigma[svd->perm[i]];
282:   MatGetSize(svd->OP,&M,&N);
283:   if (M<N) { w = u; u = v; v = w; }
284:   if (u) {
285:     SVDComputeVectors(svd);
286:     BVCopyVec(svd->U,svd->perm[i],u);
287:   }
288:   if (v) {
289:     BVCopyVec(svd->V,svd->perm[i],v);
290:   }
291:   return(0);
292: }

296: /*
297:    SVDComputeResidualNorms_Private - Computes the norms of the left and
298:    right residuals associated with the i-th computed singular triplet.
299: @*/
300: static PetscErrorCode SVDComputeResidualNorms_Private(SVD svd,PetscInt i,PetscReal *norm1,PetscReal *norm2)
301: {
303:   Vec            u,v,x = NULL,y = NULL;
304:   PetscReal      sigma;
305:   PetscInt       M,N;

308:   MatCreateVecs(svd->OP,&v,&u);
309:   SVDGetSingularTriplet(svd,i,&sigma,u,v);
310:   /* norm1 = ||A*v-sigma*u||_2 */
311:   if (norm1) {
312:     VecDuplicate(u,&x);
313:     MatMult(svd->OP,v,x);
314:     VecAXPY(x,-sigma,u);
315:     VecNorm(x,NORM_2,norm1);
316:   }
317:   /* norm2 = ||A^T*u-sigma*v||_2 */
318:   if (norm2) {
319:     VecDuplicate(v,&y);
320:     if (svd->A && svd->AT) {
321:       MatGetSize(svd->OP,&M,&N);
322:       if (M<N) {
323:         MatMult(svd->A,u,y);
324:       } else {
325:         MatMult(svd->AT,u,y);
326:       }
327:     } else {
328: #if defined(PETSC_USE_COMPLEX)
329:       MatMultHermitianTranspose(svd->OP,u,y);
330: #else
331:       MatMultTranspose(svd->OP,u,y);
332: #endif
333:     }
334:     VecAXPY(y,-sigma,v);
335:     VecNorm(y,NORM_2,norm2);
336:   }

338:   VecDestroy(&v);
339:   VecDestroy(&u);
340:   VecDestroy(&x);
341:   VecDestroy(&y);
342:   return(0);
343: }

347: /*@
348:    SVDComputeError - Computes the error (based on the residual norm) associated
349:    with the i-th singular triplet.

351:    Collective on SVD

353:    Input Parameter:
354: +  svd  - the singular value solver context
355: .  i    - the solution index
356: -  type - the type of error to compute

358:    Output Parameter:
359: .  error - the error

361:    Notes:
362:    The error can be computed in various ways, all of them based on the residual
363:    norm obtained as sqrt(n1^2+n2^2) with n1 = ||A*v-sigma*u||_2 and
364:    n2 = ||A^T*u-sigma*v||_2, where sigma is the singular value, u is the left
365:    singular vector and v is the right singular vector.

367:    Level: beginner

369: .seealso: SVDErrorType, SVDSolve()
370: @*/
371: PetscErrorCode SVDComputeError(SVD svd,PetscInt i,SVDErrorType type,PetscReal *error)
372: {
374:   PetscReal      sigma,norm1,norm2;

381:   SVDCheckSolved(svd,1);
382:   SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
383:   SVDComputeResidualNorms_Private(svd,i,&norm1,&norm2);
384:   *error = PetscSqrtReal(norm1*norm1+norm2*norm2);
385:   switch (type) {
386:     case SVD_ERROR_ABSOLUTE:
387:       break;
388:     case SVD_ERROR_RELATIVE:
389:       *error /= sigma;
390:       break;
391:     default:
392:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
393:   }
394:   return(0);
395: }