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-2009, Universidad 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

 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:   PetscTruth     flg;
 47:   PetscInt       i,*workperm;


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

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

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

 77:   PetscOptionsHasName(((PetscObject)svd)->prefix,"-svd_view",&flg);
 78:   if (flg && !PetscPreLoadingOn) { SVDView(svd,PETSC_VIEWER_STDOUT_WORLD); }

 80:   return(0);
 81: }

 85: /*@
 86:    SVDGetIterationNumber - Gets the current iteration number. If the 
 87:    call to SVDSolve() is complete, then it returns the number of iterations 
 88:    carried out by the solution method.
 89:  
 90:    Not Collective

 92:    Input Parameter:
 93: .  svd - the singular value solver context

 95:    Output Parameter:
 96: .  its - number of iterations

 98:    Level: intermediate

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

107: @*/
108: PetscErrorCode SVDGetIterationNumber(SVD svd,PetscInt *its)
109: {
113:   *its = svd->its;
114:   return(0);
115: }

119: /*@C
120:    SVDGetConvergedReason - Gets the reason why the SVDSolve() iteration was 
121:    stopped.

123:    Not Collective

125:    Input Parameter:
126: .  svd - the singular value solver context

128:    Output Parameter:
129: .  reason - negative value indicates diverged, positive value converged
130:    (see SVDConvergedReason)

132:    Possible values for reason:
133: +  SVD_CONVERGED_TOL - converged up to tolerance
134: .  SVD_DIVERGED_ITS - required more than its to reach convergence
135: -  SVD_DIVERGED_BREAKDOWN - generic breakdown in method

137:    Level: intermediate

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

141: .seealso: SVDSetTolerances(), SVDSolve(), SVDConvergedReason
142: @*/
143: PetscErrorCode SVDGetConvergedReason(SVD svd,SVDConvergedReason *reason)
144: {
148:   *reason = svd->reason;
149:   return(0);
150: }

154: /*@
155:    SVDGetConverged - Gets the number of converged singular values.

157:    Not Collective

159:    Input Parameter:
160: .  svd - the singular value solver context
161:   
162:    Output Parameter:
163: .  nconv - number of converged singular values 

165:    Note:
166:    This function should be called after SVDSolve() has finished.

168:    Level: beginner

170: @*/
171: PetscErrorCode SVDGetConverged(SVD svd,PetscInt *nconv)
172: {
176:   if (svd->reason == SVD_CONVERGED_ITERATING) {
177:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "SVDSolve must be called first");
178:   }
179:   *nconv = svd->nconv;
180:   return(0);
181: }

185: /*@
186:    SVDGetSingularTriplet - Gets the i-th triplet of the singular value decomposition
187:    as computed by SVDSolve(). The solution consists in the singular value and its left 
188:    and right singular vectors.

190:    Not Collective

192:    Input Parameters:
193: +  svd - singular value solver context 
194: -  i   - index of the solution

196:    Output Parameters:
197: +  sigma - singular value
198: .  u     - left singular vector
199: -  v     - right singular vector

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

205:    Level: beginner

207: .seealso: SVDSolve(),  SVDGetConverged()
208: @*/
209: PetscErrorCode SVDGetSingularTriplet(SVD svd, PetscInt i, PetscReal *sigma, Vec u, Vec v)
210: {
212:   PetscReal      norm;
213:   PetscInt       j,nloc,M,N;
214:   PetscScalar    *pU;
215:   Vec            w;

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

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

259:    Collective on SVD

261:    Input Parameters:
262: +  svd  - the singular value solver context
263: -  i    - the solution index

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

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

274:    Level: beginner

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

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

318:   VecDestroy(v);
319:   VecDestroy(u);
320:   if (x) { VecDestroy(x); }
321:   if (y) { VecDestroy(y); }
322:   return(0);
323: }

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

331:    Collective on SVD

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

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

343:    Level: beginner

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

355:   SVDGetSingularTriplet(svd,i,&sigma,PETSC_NULL,PETSC_NULL);
356:   SVDComputeResidualNorms(svd,i,&norm1,&norm2);
357:   *error = sqrt(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:     IPGetOperationCounters(svd->ip,dots);
392:   }
393:   return(0);
394: }