Actual source code: svdsolve.c

slepc-3.5.2 2014-10-10
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-2014, 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: /*@
 29:    SVDSolve - Solves the singular value problem.

 31:    Collective on SVD

 33:    Input Parameter:
 34: .  svd - singular value solver context obtained from SVDCreate()

 36:    Options Database Keys:
 37: +  -svd_view - print information about the solver used
 38: -  -svd_view_mat binary - save the matrix to the default binary viewer

 40:    Level: beginner

 42: .seealso: SVDCreate(), SVDSetUp(), SVDDestroy()
 43: @*/
 44: PetscErrorCode SVDSolve(SVD svd)
 45: {
 46:   PetscErrorCode    ierr;
 47:   PetscBool         flg;
 48:   PetscInt          i,*workperm;
 49:   PetscViewer       viewer;
 50:   PetscViewerFormat format;

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

 56:   /* call setup */
 57:   SVDSetUp(svd);
 58:   svd->its = 0;
 59:   svd->nconv = 0;
 60:   for (i=0;i<svd->ncv;i++) {
 61:     svd->sigma[i]  = 0.0;
 62:     svd->errest[i] = 0.0;
 63:   }
 64:   SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,svd->ncv);

 66:   (*svd->ops->solve)(svd);

 68:   /* sort singular triplets */
 69:   if (svd->which == SVD_SMALLEST) {
 70:     for (i=0;i<svd->nconv;i++) svd->perm[i] = i;
 71:     PetscSortRealWithPermutation(svd->nconv,svd->sigma,svd->perm);
 72:   } else {
 73:     PetscMalloc(sizeof(PetscInt)*svd->nconv,&workperm);
 74:     for (i=0;i<svd->nconv;i++) workperm[i] = i;
 75:     PetscSortRealWithPermutation(svd->nconv,svd->sigma,workperm);
 76:     for (i=0;i<svd->nconv;i++) svd->perm[i] = workperm[svd->nconv-i-1];
 77:     PetscFree(workperm);
 78:   }

 80:   svd->lvecsavail = (svd->leftbasis)? PETSC_TRUE: PETSC_FALSE;
 81:   PetscLogEventEnd(SVD_Solve,svd,0,0,0);

 83:   /* various viewers */
 84:   MatViewFromOptions(svd->OP,((PetscObject)svd)->prefix,"-svd_view_mat");

 86:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->prefix,"-svd_view",&viewer,&format,&flg);
 87:   if (flg && !PetscPreLoadingOn) {
 88:     PetscViewerPushFormat(viewer,format);
 89:     SVDView(svd,viewer);
 90:     PetscViewerPopFormat(viewer);
 91:     PetscViewerDestroy(&viewer);
 92:   }

 94:   /* Remove the initial subspaces */
 95:   svd->nini = 0;
 96:   svd->ninil = 0;
 97:   return(0);
 98: }

102: /*@
103:    SVDGetIterationNumber - Gets the current iteration number. If the
104:    call to SVDSolve() is complete, then it returns the number of iterations
105:    carried out by the solution method.

107:    Not Collective

109:    Input Parameter:
110: .  svd - the singular value solver context

112:    Output Parameter:
113: .  its - number of iterations

115:    Level: intermediate

117:    Notes:
118:       During the i-th iteration this call returns i-1. If SVDSolve() is
119:       complete, then parameter "its" contains either the iteration number at
120:       which convergence was successfully reached, or failure was detected.
121:       Call SVDGetConvergedReason() to determine if the solver converged or
122:       failed and why.

124: @*/
125: PetscErrorCode SVDGetIterationNumber(SVD svd,PetscInt *its)
126: {
130:   *its = svd->its;
131:   return(0);
132: }

136: /*@C
137:    SVDGetConvergedReason - Gets the reason why the SVDSolve() iteration was
138:    stopped.

140:    Not Collective

142:    Input Parameter:
143: .  svd - the singular value solver context

145:    Output Parameter:
146: .  reason - negative value indicates diverged, positive value converged
147:    (see SVDConvergedReason)

149:    Possible values for reason:
150: +  SVD_CONVERGED_TOL - converged up to tolerance
151: .  SVD_DIVERGED_ITS - required more than its to reach convergence
152: -  SVD_DIVERGED_BREAKDOWN - generic breakdown in method

154:    Level: intermediate

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

158: .seealso: SVDSetTolerances(), SVDSolve(), SVDConvergedReason
159: @*/
160: PetscErrorCode SVDGetConvergedReason(SVD svd,SVDConvergedReason *reason)
161: {
165:   *reason = svd->reason;
166:   return(0);
167: }

171: /*@
172:    SVDGetConverged - Gets the number of converged singular values.

174:    Not Collective

176:    Input Parameter:
177: .  svd - the singular value solver context

179:    Output Parameter:
180: .  nconv - number of converged singular values

182:    Note:
183:    This function should be called after SVDSolve() has finished.

185:    Level: beginner

187: @*/
188: PetscErrorCode SVDGetConverged(SVD svd,PetscInt *nconv)
189: {
193:   *nconv = svd->nconv;
194:   return(0);
195: }

199: /*@
200:    SVDGetSingularTriplet - Gets the i-th triplet of the singular value decomposition
201:    as computed by SVDSolve(). The solution consists in the singular value and its left
202:    and right singular vectors.

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

206:    Input Parameters:
207: +  svd - singular value solver context
208: -  i   - index of the solution

210:    Output Parameters:
211: +  sigma - singular value
212: .  u     - left singular vector
213: -  v     - right singular vector

215:    Note:
216:    The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
217:    Both U or V can be NULL if singular vectors are not required.

219:    Level: beginner

221: .seealso: SVDSolve(),  SVDGetConverged()
222: @*/
223: PetscErrorCode SVDGetSingularTriplet(SVD svd,PetscInt i,PetscReal *sigma,Vec u,Vec v)
224: {
226:   PetscReal      norm;
227:   PetscInt       j,M,N;
228:   Vec            w,tl,vj,uj;

234:   if (svd->reason == SVD_CONVERGED_ITERATING) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONGSTATE,"SVDSolve must be called first");
235:   if (i<0 || i>=svd->nconv) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
236:   *sigma = svd->sigma[svd->perm[i]];
237:   MatGetSize(svd->OP,&M,&N);
238:   if (M<N) { w = u; u = v; v = w; }
239:   if (u) {
240:     if (!svd->lvecsavail) {  /* generate left singular vectors on U */
241:       if (!svd->U) { SVDGetBV(svd,NULL,&svd->U); }
242:       SVDMatGetVecs(svd,NULL,&tl);
243:       BVSetSizesFromVec(svd->U,tl,svd->ncv);
244:       VecDestroy(&tl);
245:       for (j=0;j<svd->nconv;j++) {
246:         BVGetColumn(svd->V,j,&vj);
247:         BVGetColumn(svd->U,j,&uj);
248:         SVDMatMult(svd,PETSC_FALSE,vj,uj);
249:         BVRestoreColumn(svd->V,j,&vj);
250:         BVRestoreColumn(svd->U,j,&uj);
251:         BVOrthogonalizeColumn(svd->U,j,NULL,&norm,NULL);
252:         BVScaleColumn(svd->U,j,1.0/norm);
253:       }
254:       svd->lvecsavail = PETSC_TRUE;
255:     }
256:     BVCopyVec(svd->U,svd->perm[i],u);
257:   }
258:   if (v) {
259:     BVCopyVec(svd->V,svd->perm[i],v);
260:   }
261:   return(0);
262: }

266: /*@
267:    SVDComputeResidualNorms - Computes the norms of the residual vectors associated with
268:    the i-th computed singular triplet.

270:    Collective on SVD

272:    Input Parameters:
273: +  svd  - the singular value solver context
274: -  i    - the solution index

276:    Output Parameters:
277: +  norm1 - the norm ||A*v-sigma*u||_2 where sigma is the
278:            singular value, u and v are the left and right singular vectors.
279: -  norm2 - the norm ||A^T*u-sigma*v||_2 with the same sigma, u and v

281:    Note:
282:    The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
283:    Both output parameters can be NULL on input if not needed.

285:    Level: beginner

287: .seealso: SVDSolve(), SVDGetConverged(), SVDComputeRelativeError()
288: @*/
289: PetscErrorCode SVDComputeResidualNorms(SVD svd,PetscInt i,PetscReal *norm1,PetscReal *norm2)
290: {
292:   Vec            u,v,x = NULL,y = NULL;
293:   PetscReal      sigma;
294:   PetscInt       M,N;

299:   if (svd->reason == SVD_CONVERGED_ITERATING) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONGSTATE,"SVDSolve must be called first");
300:   if (i<0 || i>=svd->nconv) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");

302:   MatGetVecs(svd->OP,&v,&u);
303:   SVDGetSingularTriplet(svd,i,&sigma,u,v);
304:   if (norm1) {
305:     VecDuplicate(u,&x);
306:     MatMult(svd->OP,v,x);
307:     VecAXPY(x,-sigma,u);
308:     VecNorm(x,NORM_2,norm1);
309:   }
310:   if (norm2) {
311:     VecDuplicate(v,&y);
312:     if (svd->A && svd->AT) {
313:       MatGetSize(svd->OP,&M,&N);
314:       if (M<N) {
315:         MatMult(svd->A,u,y);
316:       } else {
317:         MatMult(svd->AT,u,y);
318:       }
319:     } else {
320: #if defined(PETSC_USE_COMPLEX)
321:       MatMultHermitianTranspose(svd->OP,u,y);
322: #else
323:       MatMultTranspose(svd->OP,u,y);
324: #endif
325:     }
326:     VecAXPY(y,-sigma,v);
327:     VecNorm(y,NORM_2,norm2);
328:   }

330:   VecDestroy(&v);
331:   VecDestroy(&u);
332:   VecDestroy(&x);
333:   VecDestroy(&y);
334:   return(0);
335: }

339: /*@
340:    SVDComputeRelativeError - Computes the relative error bound associated
341:    with the i-th singular triplet.

343:    Collective on SVD

345:    Input Parameter:
346: +  svd - the singular value solver context
347: -  i   - the solution index

349:    Output Parameter:
350: .  error - the relative error bound, computed as sqrt(n1^2+n2^2)/sigma
351:    where n1 = ||A*v-sigma*u||_2 , n2 = ||A^T*u-sigma*v||_2 , sigma is the singular value,
352:    u and v are the left and right singular vectors.
353:    If sigma is too small the relative error is computed as sqrt(n1^2+n2^2).

355:    Level: beginner

357: .seealso: SVDSolve(), SVDComputeResidualNorms()
358: @*/
359: PetscErrorCode SVDComputeRelativeError(SVD svd,PetscInt i,PetscReal *error)
360: {
362:   PetscReal      sigma,norm1,norm2;

368:   SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
369:   SVDComputeResidualNorms(svd,i,&norm1,&norm2);
370:   *error = PetscSqrtReal(norm1*norm1+norm2*norm2);
371:   if (sigma>*error) *error /= sigma;
372:   return(0);
373: }