Actual source code: svdsetup.c

slepc-3.15.1 2021-05-28
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SVD routines for setting up the solver
 12: */

 14: #include <slepc/private/svdimpl.h>

 16: /*@
 17:    SVDSetOperators - Set the matrices associated with the singular value problem.

 19:    Collective on svd

 21:    Input Parameters:
 22: +  svd - the singular value solver context
 23: .  A   - the matrix associated with the singular value problem
 24: -  B   - the second matrix in the case of GSVD

 26:    Level: beginner

 28: .seealso: SVDSolve(), SVDGetOperators()
 29: @*/
 30: PetscErrorCode SVDSetOperators(SVD svd,Mat A,Mat B)
 31: {
 33:   PetscInt       Ma,Na,Mb,Nb,ma,na,mb,nb,M0,N0,m0,n0;
 34:   PetscBool      samesize=PETSC_TRUE;


 43:   /* Check matrix sizes */
 44:   MatGetSize(A,&Ma,&Na);
 45:   MatGetLocalSize(A,&ma,&na);
 46:   if (svd->OP) {
 47:     MatGetSize(svd->OP,&M0,&N0);
 48:     MatGetLocalSize(svd->OP,&m0,&n0);
 49:     if (M0!=Ma || N0!=Na || m0!=ma || n0!=na) samesize = PETSC_FALSE;
 50:   }
 51:   if (B) {
 52:     MatGetSize(B,&Mb,&Nb);
 53:     MatGetLocalSize(B,&mb,&nb);
 54:     if (Na!=Nb) SETERRQ2(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONG,"Different number of columns in A (%D) and B (%D)",Na,Nb);
 55:     if (na!=nb) SETERRQ2(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONG,"Different local column size in A (%D) and B (%D)",na,nb);
 56:     if (svd->OPb) {
 57:       MatGetSize(svd->OPb,&M0,&N0);
 58:       MatGetLocalSize(svd->OPb,&m0,&n0);
 59:       if (M0!=Mb || N0!=Nb || m0!=mb || n0!=nb) samesize = PETSC_FALSE;
 60:     }
 61:   }

 63:   PetscObjectReference((PetscObject)A);
 64:   if (B) { PetscObjectReference((PetscObject)B); }
 65:   if (svd->state && !samesize) {
 66:     SVDReset(svd);
 67:   } else {
 68:     MatDestroy(&svd->OP);
 69:     MatDestroy(&svd->OPb);
 70:     MatDestroy(&svd->A);
 71:     MatDestroy(&svd->B);
 72:     MatDestroy(&svd->AT);
 73:     MatDestroy(&svd->BT);
 74:   }
 75:   svd->OP  = A;
 76:   svd->OPb = B;
 77:   svd->state = SVD_STATE_INITIAL;
 78:   return(0);
 79: }

 81: /*@
 82:    SVDGetOperators - Get the matrices associated with the singular value problem.

 84:    Collective on svd

 86:    Input Parameter:
 87: .  svd - the singular value solver context

 89:    Output Parameters:
 90: +  A  - the matrix associated with the singular value problem
 91: -  B  - the second matrix in the case of GSVD

 93:    Level: intermediate

 95: .seealso: SVDSolve(), SVDSetOperators()
 96: @*/
 97: PetscErrorCode SVDGetOperators(SVD svd,Mat *A,Mat *B)
 98: {
101:   if (A) *A = svd->OP;
102:   if (B) *B = svd->OPb;
103:   return(0);
104: }

106: /*@
107:    SVDSetUp - Sets up all the internal data structures necessary for the
108:    execution of the singular value solver.

110:    Collective on svd

112:    Input Parameter:
113: .  svd   - singular value solver context

115:    Notes:
116:    This function need not be called explicitly in most cases, since SVDSolve()
117:    calls it. It can be useful when one wants to measure the set-up time
118:    separately from the solve time.

120:    Level: developer

122: .seealso: SVDCreate(), SVDSolve(), SVDDestroy()
123: @*/
124: PetscErrorCode SVDSetUp(SVD svd)
125: {
127:   PetscBool      flg;
128:   PetscInt       M,N,P=0,k,maxnsol;
129:   SlepcSC        sc;
130:   Vec            *T;
131:   BV             bv;

135:   if (svd->state) return(0);
136:   PetscLogEventBegin(SVD_SetUp,svd,0,0,0);

138:   /* reset the convergence flag from the previous solves */
139:   svd->reason = SVD_CONVERGED_ITERATING;

141:   /* Set default solver type (SVDSetFromOptions was not called) */
142:   if (!((PetscObject)svd)->type_name) {
143:     SVDSetType(svd,SVDCROSS);
144:   }
145:   if (!svd->ds) { SVDGetDS(svd,&svd->ds); }

147:   /* check matrices */
148:   if (!svd->OP) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONGSTATE,"SVDSetOperators() must be called first");

150:   /* Set default problem type */
151:   if (!svd->problem_type) {
152:     if (svd->OPb) {
153:       SVDSetProblemType(svd,SVD_GENERALIZED);
154:     } else {
155:       SVDSetProblemType(svd,SVD_STANDARD);
156:     }
157:   } else if (!svd->OPb && svd->isgeneralized) {
158:     PetscInfo(svd,"Problem type set as generalized but no matrix B was provided; reverting to a standard singular value problem\n");
159:     svd->isgeneralized = PETSC_FALSE;
160:     svd->problem_type = SVD_STANDARD;
161:   } else if (svd->OPb && !svd->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_INCOMP,"Inconsistent SVD state: the problem type does not match the number of matrices");

163:   /* determine how to handle the transpose */
164:   svd->expltrans = PETSC_TRUE;
165:   if (svd->impltrans) svd->expltrans = PETSC_FALSE;
166:   else {
167:     MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);
168:     if (!flg) svd->expltrans = PETSC_FALSE;
169:     else {
170:       PetscObjectTypeCompareAny((PetscObject)svd,&flg,SVDLAPACK,SVDSCALAPACK,SVDELEMENTAL,"");
171:       if (flg) svd->expltrans = PETSC_FALSE;
172:     }
173:   }

175:   /* get matrix dimensions */
176:   MatGetSize(svd->OP,&M,&N);
177:   if (svd->isgeneralized) {
178:     MatGetSize(svd->OPb,&P,NULL);
179:     if (M+P<N) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"The case when [A;B] has less rows than columns is not supported");
180:   }

182:   /* build transpose matrix */
183:   MatDestroy(&svd->A);
184:   MatDestroy(&svd->AT);
185:   PetscObjectReference((PetscObject)svd->OP);
186:   if (svd->expltrans) {
187:     if (svd->isgeneralized || M>=N) {
188:       svd->A = svd->OP;
189:       MatHermitianTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->AT);
190:     } else {
191:       MatHermitianTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->A);
192:       svd->AT = svd->OP;
193:     }
194:   } else {
195:     if (svd->isgeneralized || M>=N) {
196:       svd->A = svd->OP;
197:       MatCreateHermitianTranspose(svd->OP,&svd->AT);
198:     } else {
199:       MatCreateHermitianTranspose(svd->OP,&svd->A);
200:       svd->AT = svd->OP;
201:     }
202:   }

204:   /* build transpose matrix B for GSVD */
205:   if (svd->isgeneralized) {
206:     MatDestroy(&svd->B);
207:     MatDestroy(&svd->BT);
208:     PetscObjectReference((PetscObject)svd->OPb);
209:     if (svd->expltrans) {
210:       svd->B = svd->OPb;
211:       MatHermitianTranspose(svd->OPb,MAT_INITIAL_MATRIX,&svd->BT);
212:     } else {
213:       svd->B = svd->OPb;
214:       MatCreateHermitianTranspose(svd->OPb,&svd->BT);
215:     }
216:   }

218:   if (!svd->isgeneralized && M<N) {
219:     /* swap initial vectors */
220:     if (svd->nini || svd->ninil) {
221:       T=svd->ISL; svd->ISL=svd->IS; svd->IS=T;
222:       k=svd->ninil; svd->ninil=svd->nini; svd->nini=k;
223:     }
224:     /* swap basis vectors */
225:     if (!svd->swapped) {  /* only the first time in case of multiple calls */
226:       bv=svd->V; svd->V=svd->U; svd->U=bv;
227:       svd->swapped = PETSC_TRUE;
228:     }
229:   }

231:   maxnsol = svd->isgeneralized? PetscMin(PetscMin(M,N),P): PetscMin(M,N);
232:   svd->ncv = PetscMin(svd->ncv,maxnsol);
233:   svd->nsv = PetscMin(svd->nsv,maxnsol);
234:   if (svd->ncv!=PETSC_DEFAULT && svd->nsv > svd->ncv) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"nsv bigger than ncv");

236:   /* call specific solver setup */
237:   (*svd->ops->setup)(svd);

239:   /* set tolerance if not yet set */
240:   if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;

242:   /* fill sorting criterion context */
243:   DSGetSlepcSC(svd->ds,&sc);
244:   sc->comparison    = (svd->which==SVD_LARGEST)? SlepcCompareLargestReal: SlepcCompareSmallestReal;
245:   sc->comparisonctx = NULL;
246:   sc->map           = NULL;
247:   sc->mapobj        = NULL;

249:   /* process initial vectors */
250:   if (svd->nini<0) {
251:     k = -svd->nini;
252:     if (k>svd->ncv) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The number of initial vectors is larger than ncv");
253:     BVInsertVecs(svd->V,0,&k,svd->IS,PETSC_TRUE);
254:     SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
255:     svd->nini = k;
256:   }
257:   if (svd->ninil<0) {
258:     k = 0;
259:     if (svd->leftbasis) {
260:       k = -svd->ninil;
261:       if (k>svd->ncv) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The number of left initial vectors is larger than ncv");
262:       BVInsertVecs(svd->U,0,&k,svd->ISL,PETSC_TRUE);
263:     } else {
264:       PetscInfo(svd,"Ignoring initial left vectors\n");
265:     }
266:     SlepcBasisDestroy_Private(&svd->ninil,&svd->ISL);
267:     svd->ninil = k;
268:   }

270:   PetscLogEventEnd(SVD_SetUp,svd,0,0,0);
271:   svd->state = SVD_STATE_SETUP;
272:   return(0);
273: }

275: /*@C
276:    SVDSetInitialSpaces - Specify two basis of vectors that constitute the initial
277:    right and/or left spaces, that is, a rough approximation to the right and/or
278:    left singular subspaces from which the solver starts to iterate.

280:    Collective on svd

282:    Input Parameter:
283: +  svd   - the singular value solver context
284: .  nr    - number of right vectors
285: .  isr   - set of basis vectors of the right initial space
286: .  nl    - number of left vectors
287: -  isl   - set of basis vectors of the left initial space

289:    Notes:
290:    It is not necessary to provide both sets of vectors.

292:    Some solvers start to iterate on a single vector (initial vector). In that case,
293:    the other vectors are ignored.

295:    These vectors do not persist from one SVDSolve() call to the other, so the
296:    initial space should be set every time.

298:    The vectors do not need to be mutually orthonormal, since they are explicitly
299:    orthonormalized internally.

301:    Common usage of this function is when the user can provide a rough approximation
302:    of the wanted singular space. Then, convergence may be faster.

304:    Level: intermediate
305: @*/
306: PetscErrorCode SVDSetInitialSpaces(SVD svd,PetscInt nr,Vec isr[],PetscInt nl,Vec isl[])
307: {

314:   if (nr<0) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument nr cannot be negative");
315:   if (nr>0) {
318:   }
319:   if (nl<0) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument nl cannot be negative");
320:   if (nl>0) {
323:   }
324:   SlepcBasisReference_Private(nr,isr,&svd->nini,&svd->IS);
325:   SlepcBasisReference_Private(nl,isl,&svd->ninil,&svd->ISL);
326:   if (nr>0 || nl>0) svd->state = SVD_STATE_INITIAL;
327:   return(0);
328: }

330: /*
331:   SVDSetDimensions_Default - Set reasonable values for ncv, mpd if not set
332:   by the user. This is called at setup.
333:  */
334: PetscErrorCode SVDSetDimensions_Default(SVD svd)
335: {
337:   PetscInt       N;

340:   MatGetSize(svd->A,NULL,&N);
341:   if (svd->ncv!=PETSC_DEFAULT) { /* ncv set */
342:     if (svd->ncv<svd->nsv) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The value of ncv must be at least nsv");
343:   } else if (svd->mpd!=PETSC_DEFAULT) { /* mpd set */
344:     svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
345:   } else { /* neither set: defaults depend on nsv being small or large */
346:     if (svd->nsv<500) svd->ncv = PetscMin(N,PetscMax(2*svd->nsv,10));
347:     else {
348:       svd->mpd = 500;
349:       svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
350:     }
351:   }
352:   if (svd->mpd==PETSC_DEFAULT) svd->mpd = svd->ncv;
353:   return(0);
354: }

356: /*@
357:    SVDAllocateSolution - Allocate memory storage for common variables such
358:    as the singular values and the basis vectors.

360:    Collective on svd

362:    Input Parameters:
363: +  svd   - eigensolver context
364: -  extra - number of additional positions, used for methods that require a
365:            working basis slightly larger than ncv

367:    Developers Notes:
368:    This is SLEPC_EXTERN because it may be required by user plugin SVD
369:    implementations.

371:    This is called at setup after setting the value of ncv and the flag leftbasis.

373:    Level: developer
374: @*/
375: PetscErrorCode SVDAllocateSolution(SVD svd,PetscInt extra)
376: {
378:   PetscInt       oldsize,requested;
379:   Vec            tr,tl;

382:   requested = svd->ncv + extra;

384:   /* oldsize is zero if this is the first time setup is called */
385:   BVGetSizes(svd->V,NULL,NULL,&oldsize);

387:   /* allocate sigma */
388:   if (requested != oldsize || !svd->sigma) {
389:     PetscFree3(svd->sigma,svd->perm,svd->errest);
390:     PetscMalloc3(requested,&svd->sigma,requested,&svd->perm,requested,&svd->errest);
391:     PetscLogObjectMemory((PetscObject)svd,PetscMax(0,requested-oldsize)*(2*sizeof(PetscReal)+sizeof(PetscInt)));
392:   }
393:   /* allocate V */
394:   if (!svd->V) { SVDGetBV(svd,&svd->V,NULL); }
395:   if (!oldsize) {
396:     if (!((PetscObject)(svd->V))->type_name) {
397:       BVSetType(svd->V,BVSVEC);
398:     }
399:     MatCreateVecsEmpty(svd->A,&tr,NULL);
400:     BVSetSizesFromVec(svd->V,tr,requested);
401:     VecDestroy(&tr);
402:   } else {
403:     BVResize(svd->V,requested,PETSC_FALSE);
404:   }
405:   /* allocate U */
406:   if (svd->leftbasis && !svd->isgeneralized) {
407:     if (!svd->U) { SVDGetBV(svd,NULL,&svd->U); }
408:     if (!oldsize) {
409:       if (!((PetscObject)(svd->U))->type_name) {
410:         BVSetType(svd->U,BVSVEC);
411:       }
412:       MatCreateVecsEmpty(svd->A,NULL,&tl);
413:       BVSetSizesFromVec(svd->U,tl,requested);
414:       VecDestroy(&tl);
415:     } else {
416:       BVResize(svd->U,requested,PETSC_FALSE);
417:     }
418:   } else if (svd->isgeneralized) {  /* left basis for the GSVD */
419:     if (!svd->U) { SVDGetBV(svd,NULL,&svd->U); }
420:     if (!oldsize) {
421:       if (!((PetscObject)(svd->U))->type_name) {
422:         BVSetType(svd->U,BVSVEC);
423:       }
424:       SVDCreateLeftTemplate(svd,&tl);
425:       BVSetSizesFromVec(svd->U,tl,requested);
426:       VecDestroy(&tl);
427:     } else {
428:       BVResize(svd->U,requested,PETSC_FALSE);
429:     }
430:   }
431:   return(0);
432: }