Actual source code: svdbasic.c
slepc-3.5.2 2014-10-10
1: /*
2: The basic SVD routines, Create, View, etc. are here.
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*/
26: PetscFunctionList SVDList = 0;
27: PetscBool SVDRegisterAllCalled = PETSC_FALSE;
28: PetscClassId SVD_CLASSID = 0;
29: PetscLogEvent SVD_SetUp = 0,SVD_Solve = 0;
33: /*@C
34: SVDView - Prints the SVD data structure.
36: Collective on SVD
38: Input Parameters:
39: + svd - the singular value solver context
40: - viewer - optional visualization context
42: Options Database Key:
43: . -svd_view - Calls SVDView() at end of SVDSolve()
45: Note:
46: The available visualization contexts include
47: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
48: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
49: output where only the first processor opens
50: the file. All other processors send their
51: data to the first processor to print.
53: The user can open an alternative visualization context with
54: PetscViewerASCIIOpen() - output to a specified file.
56: Level: beginner
58: .seealso: STView(), PetscViewerASCIIOpen()
59: @*/
60: PetscErrorCode SVDView(SVD svd,PetscViewer viewer)
61: {
63: PetscBool isascii,isshell;
67: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
71: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
72: if (isascii) {
73: PetscObjectPrintClassNamePrefixType((PetscObject)svd,viewer);
74: if (svd->ops->view) {
75: PetscViewerASCIIPushTab(viewer);
76: (*svd->ops->view)(svd,viewer);
77: PetscViewerASCIIPopTab(viewer);
78: }
79: PetscViewerASCIIPrintf(viewer," transpose mode: %s\n",svd->impltrans?"implicit":"explicit");
80: if (svd->which == SVD_LARGEST) {
81: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: largest\n");
82: } else {
83: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: smallest\n");
84: }
85: PetscViewerASCIIPrintf(viewer," number of singular values (nsv): %D\n",svd->nsv);
86: PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %D\n",svd->ncv);
87: PetscViewerASCIIPrintf(viewer," maximum dimension of projected problem (mpd): %D\n",svd->mpd);
88: PetscViewerASCIIPrintf(viewer," maximum number of iterations: %D\n",svd->max_it);
89: PetscViewerASCIIPrintf(viewer," tolerance: %g\n",(double)svd->tol);
90: if (svd->nini) {
91: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial space: %D\n",PetscAbs(svd->nini));
92: }
93: if (svd->ninil) {
94: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial left space: %D\n",PetscAbs(svd->ninil));
95: }
96: } else {
97: if (svd->ops->view) {
98: (*svd->ops->view)(svd,viewer);
99: }
100: }
101: PetscObjectTypeCompareAny((PetscObject)svd,&isshell,SVDCROSS,SVDCYCLIC,"");
102: if (!isshell) {
103: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
104: if (!svd->V) { SVDGetBV(svd,&svd->V,NULL); }
105: BVView(svd->V,viewer);
106: if (!svd->ds) { SVDGetDS(svd,&svd->ds); }
107: DSView(svd->ds,viewer);
108: PetscViewerPopFormat(viewer);
109: }
110: return(0);
111: }
115: /*@
116: SVDPrintSolution - Prints the computed singular values.
118: Collective on SVD
120: Input Parameters:
121: + svd - the singular value solver context
122: - viewer - optional visualization context
124: Options Database Key:
125: . -svd_terse - print only minimal information
127: Note:
128: By default, this function prints a table with singular values and associated
129: relative errors. With -svd_terse only the singular values are printed.
131: Level: intermediate
133: .seealso: PetscViewerASCIIOpen()
134: @*/
135: PetscErrorCode SVDPrintSolution(SVD svd,PetscViewer viewer)
136: {
137: PetscBool terse,errok,isascii;
138: PetscReal error,sigma;
139: PetscInt i,j;
144: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
147: if (!svd->sigma) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONGSTATE,"SVDSolve must be called first");
148: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
149: if (!isascii) return(0);
151: PetscOptionsHasName(NULL,"-svd_terse",&terse);
152: if (terse) {
153: if (svd->nconv<svd->nsv) {
154: PetscViewerASCIIPrintf(viewer," Problem: less than %D singular values converged\n\n",svd->nsv);
155: } else {
156: errok = PETSC_TRUE;
157: for (i=0;i<svd->nsv;i++) {
158: SVDComputeRelativeError(svd,i,&error);
159: errok = (errok && error<5.0*svd->tol)? PETSC_TRUE: PETSC_FALSE;
160: }
161: if (errok) {
162: PetscViewerASCIIPrintf(viewer," All requested singular values computed up to the required tolerance:");
163: for (i=0;i<=(svd->nsv-1)/8;i++) {
164: PetscViewerASCIIPrintf(viewer,"\n ");
165: for (j=0;j<PetscMin(8,svd->nsv-8*i);j++) {
166: SVDGetSingularTriplet(svd,8*i+j,&sigma,NULL,NULL);
167: PetscViewerASCIIPrintf(viewer,"%.5f",(double)sigma);
168: if (8*i+j+1<svd->nsv) { PetscViewerASCIIPrintf(viewer,", "); }
169: }
170: }
171: PetscViewerASCIIPrintf(viewer,"\n\n");
172: } else {
173: PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",svd->nsv);
174: }
175: }
176: } else {
177: PetscViewerASCIIPrintf(viewer," Number of converged approximate singular triplets: %D\n\n",svd->nconv);
178: if (svd->nconv>0) {
179: PetscViewerASCIIPrintf(viewer,
180: " sigma relative error\n"
181: " --------------------- ------------------\n");
182: for (i=0;i<svd->nconv;i++) {
183: SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
184: SVDComputeRelativeError(svd,i,&error);
185: PetscViewerASCIIPrintf(viewer," % 6f %12g\n",(double)sigma,(double)error);
186: }
187: PetscViewerASCIIPrintf(viewer,"\n");
188: }
189: }
190: return(0);
191: }
195: /*@C
196: SVDCreate - Creates the default SVD context.
198: Collective on MPI_Comm
200: Input Parameter:
201: . comm - MPI communicator
203: Output Parameter:
204: . svd - location to put the SVD context
206: Note:
207: The default SVD type is SVDCROSS
209: Level: beginner
211: .seealso: SVDSetUp(), SVDSolve(), SVDDestroy(), SVD
212: @*/
213: PetscErrorCode SVDCreate(MPI_Comm comm,SVD *outsvd)
214: {
216: SVD svd;
220: *outsvd = 0;
221: SVDInitializePackage();
222: SlepcHeaderCreate(svd,_p_SVD,struct _SVDOps,SVD_CLASSID,"SVD","Singular Value Decomposition","SVD",comm,SVDDestroy,SVDView);
224: svd->OP = NULL;
225: svd->max_it = 0;
226: svd->nsv = 1;
227: svd->ncv = 0;
228: svd->mpd = 0;
229: svd->nini = 0;
230: svd->ninil = 0;
231: svd->tol = PETSC_DEFAULT;
232: svd->which = SVD_LARGEST;
233: svd->impltrans = PETSC_FALSE;
234: svd->trackall = PETSC_FALSE;
236: svd->numbermonitors = 0;
238: svd->ds = NULL;
239: svd->U = NULL;
240: svd->V = NULL;
241: svd->rand = NULL;
242: svd->A = NULL;
243: svd->AT = NULL;
244: svd->IS = NULL;
245: svd->ISL = NULL;
246: svd->sigma = NULL;
247: svd->perm = NULL;
248: svd->errest = NULL;
249: svd->data = NULL;
251: svd->nconv = 0;
252: svd->its = 0;
253: svd->leftbasis = PETSC_FALSE;
254: svd->lvecsavail = PETSC_FALSE;
255: svd->setupcalled = 0;
256: svd->reason = SVD_CONVERGED_ITERATING;
258: PetscNewLog(svd,&svd->sc);
259: PetscRandomCreate(comm,&svd->rand);
260: PetscRandomSetSeed(svd->rand,0x12345678);
261: PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->rand);
262: *outsvd = svd;
263: return(0);
264: }
268: /*@
269: SVDReset - Resets the SVD context to the setupcalled=0 state and removes any
270: allocated objects.
272: Collective on SVD
274: Input Parameter:
275: . svd - singular value solver context obtained from SVDCreate()
277: Level: advanced
279: .seealso: SVDDestroy()
280: @*/
281: PetscErrorCode SVDReset(SVD svd)
282: {
284: PetscInt ncols;
288: if (svd->ops->reset) { (svd->ops->reset)(svd); }
289: if (svd->ds) { DSReset(svd->ds); }
290: MatDestroy(&svd->OP);
291: MatDestroy(&svd->A);
292: MatDestroy(&svd->AT);
293: BVGetSizes(svd->V,NULL,NULL,&ncols);
294: if (ncols) {
295: PetscFree3(svd->sigma,svd->perm,svd->errest);
296: }
297: BVDestroy(&svd->U);
298: BVDestroy(&svd->V);
299: svd->setupcalled = 0;
300: return(0);
301: }
305: /*@C
306: SVDDestroy - Destroys the SVD context.
308: Collective on SVD
310: Input Parameter:
311: . svd - singular value solver context obtained from SVDCreate()
313: Level: beginner
315: .seealso: SVDCreate(), SVDSetUp(), SVDSolve()
316: @*/
317: PetscErrorCode SVDDestroy(SVD *svd)
318: {
322: if (!*svd) return(0);
324: if (--((PetscObject)(*svd))->refct > 0) { *svd = 0; return(0); }
325: SVDReset(*svd);
326: if ((*svd)->ops->destroy) { (*(*svd)->ops->destroy)(*svd); }
327: DSDestroy(&(*svd)->ds);
328: PetscRandomDestroy(&(*svd)->rand);
329: PetscFree((*svd)->sc);
330: /* just in case the initial vectors have not been used */
331: SlepcBasisDestroy_Private(&(*svd)->nini,&(*svd)->IS);
332: SlepcBasisDestroy_Private(&(*svd)->ninil,&(*svd)->ISL);
333: SVDMonitorCancel(*svd);
334: PetscHeaderDestroy(svd);
335: return(0);
336: }
340: /*@C
341: SVDSetType - Selects the particular solver to be used in the SVD object.
343: Logically Collective on SVD
345: Input Parameters:
346: + svd - the singular value solver context
347: - type - a known method
349: Options Database Key:
350: . -svd_type <method> - Sets the method; use -help for a list
351: of available methods
353: Notes:
354: See "slepc/include/slepcsvd.h" for available methods. The default
355: is SVDCROSS.
357: Normally, it is best to use the SVDSetFromOptions() command and
358: then set the SVD type from the options database rather than by using
359: this routine. Using the options database provides the user with
360: maximum flexibility in evaluating the different available methods.
361: The SVDSetType() routine is provided for those situations where it
362: is necessary to set the iterative solver independently of the command
363: line or options database.
365: Level: intermediate
367: .seealso: SVDType
368: @*/
369: PetscErrorCode SVDSetType(SVD svd,SVDType type)
370: {
371: PetscErrorCode ierr,(*r)(SVD);
372: PetscBool match;
378: PetscObjectTypeCompare((PetscObject)svd,type,&match);
379: if (match) return(0);
381: PetscFunctionListFind(SVDList,type,&r);
382: if (!r) SETERRQ1(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown SVD type given: %s",type);
384: if (svd->ops->destroy) { (*svd->ops->destroy)(svd); }
385: PetscMemzero(svd->ops,sizeof(struct _SVDOps));
387: svd->setupcalled = 0;
388: PetscObjectChangeTypeName((PetscObject)svd,type);
389: (*r)(svd);
390: return(0);
391: }
395: /*@C
396: SVDGetType - Gets the SVD type as a string from the SVD object.
398: Not Collective
400: Input Parameter:
401: . svd - the singular value solver context
403: Output Parameter:
404: . name - name of SVD method
406: Level: intermediate
408: .seealso: SVDSetType()
409: @*/
410: PetscErrorCode SVDGetType(SVD svd,SVDType *type)
411: {
415: *type = ((PetscObject)svd)->type_name;
416: return(0);
417: }
421: /*@C
422: SVDRegister - Adds a method to the singular value solver package.
424: Not Collective
426: Input Parameters:
427: + name - name of a new user-defined solver
428: - function - routine to create the solver context
430: Notes:
431: SVDRegister() may be called multiple times to add several user-defined solvers.
433: Sample usage:
434: .vb
435: SVDRegister("my_solver",MySolverCreate);
436: .ve
438: Then, your solver can be chosen with the procedural interface via
439: $ SVDSetType(svd,"my_solver")
440: or at runtime via the option
441: $ -svd_type my_solver
443: Level: advanced
445: .seealso: SVDRegisterAll()
446: @*/
447: PetscErrorCode SVDRegister(const char *name,PetscErrorCode (*function)(SVD))
448: {
452: PetscFunctionListAdd(&SVDList,name,function);
453: return(0);
454: }
458: /*@
459: SVDSetBV - Associates basis vectors objects to the singular value solver.
461: Collective on SVD
463: Input Parameters:
464: + svd - singular value solver context obtained from SVDCreate()
465: . V - the basis vectors object for right singular vectors
466: - U - the basis vectors object for left singular vectors
468: Note:
469: Use SVDGetBV() to retrieve the basis vectors contexts (for example,
470: to free them at the end of the computations).
472: Level: advanced
474: .seealso: SVDGetBV()
475: @*/
476: PetscErrorCode SVDSetBV(SVD svd,BV V,BV U)
477: {
482: if (V) {
485: PetscObjectReference((PetscObject)V);
486: BVDestroy(&svd->V);
487: svd->V = V;
488: PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->V);
489: }
490: if (U) {
493: PetscObjectReference((PetscObject)U);
494: BVDestroy(&svd->U);
495: svd->U = U;
496: PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->U);
497: }
498: return(0);
499: }
503: /*@C
504: SVDGetBV - Obtain the basis vectors objects associated to the singular
505: value solver object.
507: Not Collective
509: Input Parameters:
510: . svd - singular value solver context obtained from SVDCreate()
512: Output Parameter:
513: + V - basis vectors context for right singular vectors
514: - U - basis vectors context for left singular vectors
516: Level: advanced
518: .seealso: SVDSetBV()
519: @*/
520: PetscErrorCode SVDGetBV(SVD svd,BV *V,BV *U)
521: {
526: if (V) {
527: if (!svd->V) {
528: BVCreate(PetscObjectComm((PetscObject)svd),&svd->V);
529: PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->V);
530: }
531: *V = svd->V;
532: }
533: if (U) {
534: if (!svd->U) {
535: BVCreate(PetscObjectComm((PetscObject)svd),&svd->U);
536: PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->U);
537: }
538: *U = svd->U;
539: }
540: return(0);
541: }
545: /*@
546: SVDSetDS - Associates a direct solver object to the singular value solver.
548: Collective on SVD
550: Input Parameters:
551: + svd - singular value solver context obtained from SVDCreate()
552: - ds - the direct solver object
554: Note:
555: Use SVDGetDS() to retrieve the direct solver context (for example,
556: to free it at the end of the computations).
558: Level: advanced
560: .seealso: SVDGetDS()
561: @*/
562: PetscErrorCode SVDSetDS(SVD svd,DS ds)
563: {
570: PetscObjectReference((PetscObject)ds);
571: DSDestroy(&svd->ds);
572: svd->ds = ds;
573: PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->ds);
574: return(0);
575: }
579: /*@C
580: SVDGetDS - Obtain the direct solver object associated to the singular value
581: solver object.
583: Not Collective
585: Input Parameters:
586: . svd - singular value solver context obtained from SVDCreate()
588: Output Parameter:
589: . ds - direct solver context
591: Level: advanced
593: .seealso: SVDSetDS()
594: @*/
595: PetscErrorCode SVDGetDS(SVD svd,DS *ds)
596: {
602: if (!svd->ds) {
603: DSCreate(PetscObjectComm((PetscObject)svd),&svd->ds);
604: PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->ds);
605: }
606: *ds = svd->ds;
607: return(0);
608: }