Actual source code: svdsolve.c

  1: /*
  2:       SVD routines related to the solution process.

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

  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.

 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 <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:
 37: .   -svd_view - print information about the solver used

 39:    Level: beginner

 41: .seealso: SVDCreate(), SVDSetUp(), SVDDestroy() 
 42: @*/
 43: PetscErrorCode SVDSolve(SVD svd)
 44: {
 46:   PetscBool      flg;
 47:   PetscInt       i,*workperm;
 48:   char           filename[PETSC_MAX_PATH_LEN];
 49:   PetscViewer    viewer;

 53:   if (!svd->setupcalled) { SVDSetUp(svd); }
 54:   svd->its = 0;
 55:   svd->nconv = 0;
 56:   svd->reason = SVD_CONVERGED_ITERATING;
 57:   for (i=0;i<svd->ncv;i++) svd->sigma[i]=svd->errest[i]=0.0;
 58:   SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,svd->ncv);

 60:   PetscLogEventBegin(SVD_Solve,svd,0,0,0);
 61:   (*svd->ops->solve)(svd);
 62:   PetscLogEventEnd(SVD_Solve,svd,0,0,0);

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

 76:   PetscOptionsGetString(((PetscObject)svd)->prefix,"-svd_view",filename,PETSC_MAX_PATH_LEN,&flg);
 77:   if (flg && !PetscPreLoadingOn) {
 78:     PetscViewerASCIIOpen(((PetscObject)svd)->comm,filename,&viewer);
 79:     SVDView(svd,viewer);
 80:     PetscViewerDestroy(&viewer);
 81:   }

 83:   /* Remove the initial subspace */
 84:   svd->nini = 0;
 85:   return(0);
 86: }

 90: /*@
 91:    SVDGetIterationNumber - Gets the current iteration number. If the 
 92:    call to SVDSolve() is complete, then it returns the number of iterations 
 93:    carried out by the solution method.
 94:  
 95:    Not Collective

 97:    Input Parameter:
 98: .  svd - the singular value solver context

100:    Output Parameter:
101: .  its - number of iterations

103:    Level: intermediate

105:    Notes:
106:       During the i-th iteration this call returns i-1. If SVDSolve() is 
107:       complete, then parameter "its" contains either the iteration number at
108:       which convergence was successfully reached, or failure was detected.  
109:       Call SVDGetConvergedReason() to determine if the solver converged or 
110:       failed and why.

112: @*/
113: PetscErrorCode SVDGetIterationNumber(SVD svd,PetscInt *its)
114: {
118:   *its = svd->its;
119:   return(0);
120: }

124: /*@C
125:    SVDGetConvergedReason - Gets the reason why the SVDSolve() iteration was 
126:    stopped.

128:    Not Collective

130:    Input Parameter:
131: .  svd - the singular value solver context

133:    Output Parameter:
134: .  reason - negative value indicates diverged, positive value converged
135:    (see SVDConvergedReason)

137:    Possible values for reason:
138: +  SVD_CONVERGED_TOL - converged up to tolerance
139: .  SVD_DIVERGED_ITS - required more than its to reach convergence
140: -  SVD_DIVERGED_BREAKDOWN - generic breakdown in method

142:    Level: intermediate

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

146: .seealso: SVDSetTolerances(), SVDSolve(), SVDConvergedReason
147: @*/
148: PetscErrorCode SVDGetConvergedReason(SVD svd,SVDConvergedReason *reason)
149: {
153:   *reason = svd->reason;
154:   return(0);
155: }

159: /*@
160:    SVDGetConverged - Gets the number of converged singular values.

162:    Not Collective

164:    Input Parameter:
165: .  svd - the singular value solver context
166:   
167:    Output Parameter:
168: .  nconv - number of converged singular values 

170:    Note:
171:    This function should be called after SVDSolve() has finished.

173:    Level: beginner

175: @*/
176: PetscErrorCode SVDGetConverged(SVD svd,PetscInt *nconv)
177: {
181:   if (svd->reason == SVD_CONVERGED_ITERATING) {
182:     SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_WRONGSTATE,"SVDSolve must be called first");
183:   }
184:   *nconv = svd->nconv;
185:   return(0);
186: }

190: /*@
191:    SVDGetSingularTriplet - Gets the i-th triplet of the singular value decomposition
192:    as computed by SVDSolve(). The solution consists in the singular value and its left 
193:    and right singular vectors.

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

197:    Input Parameters:
198: +  svd - singular value solver context 
199: -  i   - index of the solution

201:    Output Parameters:
202: +  sigma - singular value
203: .  u     - left singular vector
204: -  v     - right singular vector

206:    Note:
207:    The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
208:    Both U or V can be PETSC_NULL if singular vectors are not required. 

210:    Level: beginner

212: .seealso: SVDSolve(),  SVDGetConverged()
213: @*/
214: PetscErrorCode SVDGetSingularTriplet(SVD svd,PetscInt i,PetscReal *sigma,Vec u,Vec v)
215: {
217:   PetscReal      norm;
218:   PetscInt       j,M,N;
219:   Vec            w;

225:   if (svd->reason == SVD_CONVERGED_ITERATING) {
226:     SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_WRONGSTATE,"SVDSolve must be called first");
227:   }
228:   if (i<0 || i>=svd->nconv) {
229:     SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
230:   }
231:   *sigma = svd->sigma[svd->perm[i]];
232:   MatGetSize(svd->OP,&M,&N);
233:   if (M<N) { w = u; u = v; v = w; }
234:   if (u) {
235:     if (!svd->U) {
236:       VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
237:       for (j=0;j<svd->nconv;j++) {
238:         SVDMatMult(svd,PETSC_FALSE,svd->V[j],svd->U[j]);
239:         IPOrthogonalize(svd->ip,0,PETSC_NULL,j,PETSC_NULL,svd->U,svd->U[j],PETSC_NULL,&norm,PETSC_NULL);
240:         VecScale(svd->U[j],1.0/norm);
241:       }
242:     }
243:     VecCopy(svd->U[svd->perm[i]],u);
244:   }
245:   if (v) {
246:     VecCopy(svd->V[svd->perm[i]],v);
247:   }
248:   return(0);
249: }

253: /*@
254:    SVDComputeResidualNorms - Computes the norms of the residual vectors associated with 
255:    the i-th computed singular triplet.

257:    Collective on SVD

259:    Input Parameters:
260: +  svd  - the singular value solver context
261: -  i    - the solution index

263:    Output Parameters:
264: +  norm1 - the norm ||A*v-sigma*u||_2 where sigma is the 
265:            singular value, u and v are the left and right singular vectors. 
266: -  norm2 - the norm ||A^T*u-sigma*v||_2 with the same sigma, u and v

268:    Note:
269:    The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
270:    Both output parameters can be PETSC_NULL on input if not needed.

272:    Level: beginner

274: .seealso: SVDSolve(), SVDGetConverged(), SVDComputeRelativeError()
275: @*/
276: PetscErrorCode SVDComputeResidualNorms(SVD svd,PetscInt i,PetscReal *norm1,PetscReal *norm2)
277: {
279:   Vec            u,v,x = PETSC_NULL,y = PETSC_NULL;
280:   PetscReal      sigma;
281:   PetscInt       M,N;

286:   if (svd->reason == SVD_CONVERGED_ITERATING) {
287:     SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_WRONGSTATE,"SVDSolve must be called first");
288:   }
289:   if (i<0 || i>=svd->nconv) {
290:     SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
291:   }
292: 
293:   MatGetVecs(svd->OP,&v,&u);
294:   SVDGetSingularTriplet(svd,i,&sigma,u,v);
295:   if (norm1) {
296:     VecDuplicate(u,&x);
297:     MatMult(svd->OP,v,x);
298:     VecAXPY(x,-sigma,u);
299:     VecNorm(x,NORM_2,norm1);
300:   }
301:   if (norm2) {
302:     VecDuplicate(v,&y);
303:     if (svd->A && svd->AT) {
304:       MatGetSize(svd->OP,&M,&N);
305:       if (M<N) {
306:         MatMult(svd->A,u,y);
307:       } else {
308:         MatMult(svd->AT,u,y);
309:       }
310:     } else {
311:       MatMultTranspose(svd->OP,u,y);
312:     }
313:     VecAXPY(y,-sigma,v);
314:     VecNorm(y,NORM_2,norm2);
315:   }

317:   VecDestroy(&v);
318:   VecDestroy(&u);
319:   VecDestroy(&x);
320:   VecDestroy(&y);
321:   return(0);
322: }

326: /*@
327:    SVDComputeRelativeError - Computes the relative error bound associated 
328:    with the i-th singular triplet.

330:    Collective on SVD

332:    Input Parameter:
333: +  svd - the singular value solver context
334: -  i   - the solution index

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

342:    Level: beginner

344: .seealso: SVDSolve(), SVDComputeResidualNorms()
345: @*/
346: PetscErrorCode SVDComputeRelativeError(SVD svd,PetscInt i,PetscReal *error)
347: {
349:   PetscReal      sigma,norm1,norm2;

355:   SVDGetSingularTriplet(svd,i,&sigma,PETSC_NULL,PETSC_NULL);
356:   SVDComputeResidualNorms(svd,i,&norm1,&norm2);
357:   *error = PetscSqrtReal(norm1*norm1+norm2*norm2);
358:   if (sigma>*error) *error /= sigma;
359:   return(0);
360: }

364: /*@
365:    SVDGetOperationCounters - Gets the total number of matrix vector and dot 
366:    products used by the SVD object during the last SVDSolve() call.

368:    Not Collective

370:    Input Parameter:
371: .  svd - SVD context

373:    Output Parameter:
374: +  matvecs - number of matrix vector product operations
375: -  dots    - number of dot product operations

377:    Notes:
378:    These counters are reset to zero at each successive call to SVDSolve().

380:    Level: intermediate

382: @*/
383: PetscErrorCode SVDGetOperationCounters(SVD svd,PetscInt* matvecs,PetscInt* dots)
384: {
386: 
389:   if (matvecs) *matvecs = svd->matvecs;
390:   if (dots) {
391:     if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
392:     IPGetOperationCounters(svd->ip,dots);
393:   }
394:   return(0);
395: }