Actual source code: svdbasic.c
1: /*
2: The basic SVD routines, Create, View, etc. are here.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2010, 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
26: PetscFList SVDList = 0;
27: PetscCookie SVD_COOKIE = 0;
28: PetscLogEvent SVD_SetUp = 0, SVD_Solve = 0, SVD_Dense = 0;
29: static PetscTruth SVDPackageInitialized = PETSC_FALSE;
33: /*@C
34: SVDFinalizePackage - This function destroys everything in the Slepc interface to the SVD package. It is
35: called from SlepcFinalize().
37: Level: developer
39: .seealso: SlepcFinalize()
40: @*/
41: PetscErrorCode SVDFinalizePackage(void)
42: {
44: SVDPackageInitialized = PETSC_FALSE;
45: SVDList = 0;
46: return(0);
47: }
51: /*@C
52: SVDInitializePackage - This function initializes everything in the SVD package. It is called
53: from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to SVDCreate()
54: when using static libraries.
56: Input Parameter:
57: path - The dynamic library path, or PETSC_NULL
59: Level: developer
61: .seealso: SlepcInitialize()
62: @*/
63: PetscErrorCode SVDInitializePackage(char *path)
64: {
65: char logList[256];
66: char *className;
67: PetscTruth opt;
68: PetscErrorCode ierr;
71: if (SVDPackageInitialized) return(0);
72: SVDPackageInitialized = PETSC_TRUE;
73: /* Register Classes */
74: PetscCookieRegister("Singular Value Solver",&SVD_COOKIE);
75: /* Register Constructors */
76: SVDRegisterAll(path);
77: /* Register Events */
78: PetscLogEventRegister("SVDSetUp",SVD_COOKIE,&SVD_SetUp);
79: PetscLogEventRegister("SVDSolve",SVD_COOKIE,&SVD_Solve);
80: PetscLogEventRegister("SVDDense",SVD_COOKIE,&SVD_Dense);
81: /* Process info exclusions */
82: PetscOptionsGetString(PETSC_NULL, "-log_info_exclude", logList, 256, &opt);
83: if (opt) {
84: PetscStrstr(logList, "svd", &className);
85: if (className) {
86: PetscInfoDeactivateClass(SVD_COOKIE);
87: }
88: }
89: /* Process summary exclusions */
90: PetscOptionsGetString(PETSC_NULL, "-log_summary_exclude", logList, 256, &opt);
91: if (opt) {
92: PetscStrstr(logList, "svd", &className);
93: if (className) {
94: PetscLogEventDeactivateClass(SVD_COOKIE);
95: }
96: }
97: PetscRegisterFinalize(SVDFinalizePackage);
98: return(0);
99: }
103: /*@C
104: SVDView - Prints the SVD data structure.
106: Collective on SVD
108: Input Parameters:
109: + svd - the singular value solver context
110: - viewer - optional visualization context
112: Options Database Key:
113: . -svd_view - Calls SVDView() at end of SVDSolve()
115: Note:
116: The available visualization contexts include
117: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
118: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
119: output where only the first processor opens
120: the file. All other processors send their
121: data to the first processor to print.
123: The user can open an alternative visualization context with
124: PetscViewerASCIIOpen() - output to a specified file.
126: Level: beginner
128: .seealso: STView(), PetscViewerASCIIOpen()
129: @*/
130: PetscErrorCode SVDView(SVD svd,PetscViewer viewer)
131: {
133: const SVDType type;
134: PetscTruth isascii;
138: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(((PetscObject)svd)->comm);
142: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
143: if (isascii) {
144: PetscViewerASCIIPrintf(viewer,"SVD Object:\n");
145: SVDGetType(svd,&type);
146: if (type) {
147: PetscViewerASCIIPrintf(viewer," method: %s\n",type);
148: } else {
149: PetscViewerASCIIPrintf(viewer," method: not yet set\n");
150: }
151: switch (svd->transmode) {
152: case SVD_TRANSPOSE_EXPLICIT:
153: PetscViewerASCIIPrintf(viewer," transpose mode: explicit\n");
154: break;
155: case SVD_TRANSPOSE_IMPLICIT:
156: PetscViewerASCIIPrintf(viewer," transpose mode: implicit\n");
157: break;
158: default:
159: PetscViewerASCIIPrintf(viewer," transpose mode: not yet set\n");
160: }
161: if (svd->which == SVD_LARGEST) {
162: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: largest\n");
163: } else {
164: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: smallest\n");
165: }
166: PetscViewerASCIIPrintf(viewer," number of singular values (nsv): %d\n",svd->nsv);
167: PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %d\n",svd->ncv);
168: PetscViewerASCIIPrintf(viewer," maximum dimension of projected problem (mpd): %d\n",svd->mpd);
169: PetscViewerASCIIPrintf(viewer," maximum number of iterations: %d\n",svd->max_it);
170: PetscViewerASCIIPrintf(viewer," tolerance: %g\n",svd->tol);
171: if (svd->nini!=0) {
172: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial space: %d\n",PetscAbs(svd->nini));
173: }
174: if (svd->ops->view) {
175: PetscViewerASCIIPushTab(viewer);
176: (*svd->ops->view)(svd,viewer);
177: PetscViewerASCIIPopTab(viewer);
178: }
179: PetscViewerASCIIPushTab(viewer);
180: IPView(svd->ip,viewer);
181: PetscViewerASCIIPopTab(viewer);
182: } else {
183: if (svd->ops->view) {
184: (*svd->ops->view)(svd,viewer);
185: }
186: }
187: return(0);
188: }
192: /*@C
193: SVDCreate - Creates the default SVD context.
195: Collective on MPI_Comm
197: Input Parameter:
198: . comm - MPI communicator
200: Output Parameter:
201: . svd - location to put the SVD context
203: Note:
204: The default SVD type is SVDCROSS
206: Level: beginner
208: .seealso: SVDSetUp(), SVDSolve(), SVDDestroy(), SVD
209: @*/
210: PetscErrorCode SVDCreate(MPI_Comm comm,SVD *outsvd)
211: {
213: SVD svd;
218: PetscHeaderCreate(svd,_p_SVD,struct _SVDOps,SVD_COOKIE,-1,"SVD",comm,SVDDestroy,SVDView);
219: *outsvd = svd;
221: PetscMemzero(svd->ops,sizeof(struct _SVDOps));
223: svd->OP = PETSC_NULL;
224: svd->A = PETSC_NULL;
225: svd->AT = PETSC_NULL;
226: svd->transmode = (SVDTransposeMode)PETSC_DECIDE;
227: svd->sigma = PETSC_NULL;
228: svd->perm = PETSC_NULL;
229: svd->U = PETSC_NULL;
230: svd->V = PETSC_NULL;
231: svd->IS = PETSC_NULL;
232: svd->which = SVD_LARGEST;
233: svd->n = 0;
234: svd->nconv = 0;
235: svd->nsv = 1;
236: svd->ncv = 0;
237: svd->mpd = 0;
238: svd->nini = 0;
239: svd->its = 0;
240: svd->max_it = 0;
241: svd->tol = 1e-7;
242: svd->errest = PETSC_NULL;
243: svd->data = PETSC_NULL;
244: svd->setupcalled = 0;
245: svd->reason = SVD_CONVERGED_ITERATING;
246: svd->numbermonitors = 0;
247: svd->matvecs = 0;
248: svd->trackall = PETSC_FALSE;
250: IPCreate(comm,&svd->ip);
251: IPSetOptionsPrefix(svd->ip,((PetscObject)svd)->prefix);
252: IPAppendOptionsPrefix(svd->ip,"svd_");
253: PetscLogObjectParent(svd,svd->ip);
255: PetscPublishAll(svd);
256: return(0);
257: }
258:
261: /*@
262: SVDDestroy - Destroys the SVD context.
264: Collective on SVD
266: Input Parameter:
267: . svd - singular value solver context obtained from SVDCreate()
269: Level: beginner
271: .seealso: SVDCreate(), SVDSetUp(), SVDSolve()
272: @*/
273: PetscErrorCode SVDDestroy(SVD svd)
274: {
276: PetscInt i;
277: PetscScalar *p;
278:
281: if (--((PetscObject)svd)->refct > 0) return(0);
283: /* if memory was published with AMS then destroy it */
284: PetscObjectDepublish(svd);
286: if (svd->ops->destroy) {
287: (*svd->ops->destroy)(svd);
288: }
290: if (svd->OP) { MatDestroy(svd->OP); }
291: if (svd->A) { MatDestroy(svd->A); }
292: if (svd->AT) { MatDestroy(svd->AT); }
293: if (svd->n) {
294: PetscFree(svd->sigma);
295: PetscFree(svd->perm);
296: PetscFree(svd->errest);
297: if (svd->U) {
298: VecGetArray(svd->U[0],&p);
299: for (i=0;i<svd->n;i++) {
300: VecDestroy(svd->U[i]);
301: }
302: PetscFree(p);
303: PetscFree(svd->U);
304: }
305: VecGetArray(svd->V[0],&p);
306: for (i=0;i<svd->n;i++) {
307: VecDestroy(svd->V[i]);
308: }
309: PetscFree(p);
310: PetscFree(svd->V);
311: }
312: SVDMonitorCancel(svd);
313:
314: IPDestroy(svd->ip);
315: if (svd->rand) {
316: PetscRandomDestroy(svd->rand);
317: }
318:
319: PetscHeaderDestroy(svd);
320: return(0);
321: }
325: PetscErrorCode SVDDestroy_Default(SVD svd)
326: {
331: PetscFree(svd->data);
332: return(0);
333: }
337: /*@C
338: SVDSetType - Selects the particular solver to be used in the SVD object.
340: Collective on SVD
342: Input Parameters:
343: + svd - the singular value solver context
344: - type - a known method
346: Options Database Key:
347: . -svd_type <method> - Sets the method; use -help for a list
348: of available methods
349:
350: Notes:
351: See "slepc/include/slepcsvd.h" for available methods. The default
352: is SVDCROSS.
354: Normally, it is best to use the SVDSetFromOptions() command and
355: then set the SVD type from the options database rather than by using
356: this routine. Using the options database provides the user with
357: maximum flexibility in evaluating the different available methods.
358: The SVDSetType() routine is provided for those situations where it
359: is necessary to set the iterative solver independently of the command
360: line or options database.
362: Level: intermediate
364: .seealso: SVDType
365: @*/
366: PetscErrorCode SVDSetType(SVD svd,const SVDType type)
367: {
368: PetscErrorCode ierr,(*r)(SVD);
369: PetscTruth match;
375: PetscTypeCompare((PetscObject)svd,type,&match);
376: if (match) return(0);
378: if (svd->data) {
379: /* destroy the old private SVD context */
380: (*svd->ops->destroy)(svd);
381: svd->data = 0;
382: }
384: PetscFListFind(SVDList,((PetscObject)svd)->comm,type,(void (**)(void)) &r);
386: if (!r) SETERRQ1(1,"Unknown SVD type given: %s",type);
388: svd->setupcalled = 0;
389: PetscMemzero(svd->ops,sizeof(struct _SVDOps));
390: (*r)(svd);
392: PetscObjectChangeTypeName((PetscObject)svd,type);
393: return(0);
394: }
398: /*@C
399: SVDGetType - Gets the SVD type as a string from the SVD object.
401: Not Collective
403: Input Parameter:
404: . svd - the singular value solver context
406: Output Parameter:
407: . name - name of SVD method
409: Level: intermediate
411: .seealso: SVDSetType()
412: @*/
413: PetscErrorCode SVDGetType(SVD svd,const SVDType *type)
414: {
418: *type = ((PetscObject)svd)->type_name;
419: return(0);
420: }
422: /*MC
423: SVDRegisterDynamic - Adds a method to the singular value solver package.
425: Synopsis:
426: SVDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(SVD))
428: Not Collective
430: Input Parameters:
431: + name_solver - name of a new user-defined solver
432: . path - path (either absolute or relative) the library containing this solver
433: . name_create - name of routine to create the solver context
434: - routine_create - routine to create the solver context
436: Notes:
437: SVDRegisterDynamic() may be called multiple times to add several user-defined solvers.
439: If dynamic libraries are used, then the fourth input argument (routine_create)
440: is ignored.
442: Sample usage:
443: .vb
444: SVDRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
445: "MySolverCreate",MySolverCreate);
446: .ve
448: Then, your solver can be chosen with the procedural interface via
449: $ SVDSetType(svd,"my_solver")
450: or at runtime via the option
451: $ -svd_type my_solver
453: Level: advanced
455: Environmental variables such as ${PETSC_ARCH}, ${SLEPC_DIR},
456: and others of the form ${any_environmental_variable} occuring in pathname will be
457: replaced with appropriate values.
459: .seealso: SVDRegisterDestroy(), SVDRegisterAll()
461: M*/
465: /*@C
466: SVDRegister - See SVDRegisterDynamic()
468: Level: advanced
469: @*/
470: PetscErrorCode SVDRegister(const char *sname,const char *path,const char *name,PetscErrorCode (*function)(SVD))
471: {
473: char fullname[256];
476: PetscFListConcat(path,name,fullname);
477: PetscFListAdd(&SVDList,sname,fullname,(void (*)(void))function);
478: return(0);
479: }
483: /*@
484: SVDRegisterDestroy - Frees the list of SVD methods that were
485: registered by SVDRegisterDynamic().
487: Not Collective
489: Level: advanced
491: .seealso: SVDRegisterDynamic(), SVDRegisterAll()
492: @*/
493: PetscErrorCode SVDRegisterDestroy(void)
494: {
498: PetscFListDestroy(&SVDList);
499: SVDRegisterAll(PETSC_NULL);
500: return(0);
501: }
505: /*@
506: SVDSetIP - Associates an inner product object to the
507: singular value solver.
509: Collective on SVD
511: Input Parameters:
512: + svd - singular value solver context obtained from SVDCreate()
513: - ip - the inner product object
515: Note:
516: Use SVDGetIP() to retrieve the inner product context (for example,
517: to free it at the end of the computations).
519: Level: advanced
521: .seealso: SVDGetIP()
522: @*/
523: PetscErrorCode SVDSetIP(SVD svd,IP ip)
524: {
531: PetscObjectReference((PetscObject)ip);
532: IPDestroy(svd->ip);
533: svd->ip = ip;
534: return(0);
535: }
539: /*@C
540: SVDGetIP - Obtain the inner product object associated
541: to the singular value solver object.
543: Not Collective
545: Input Parameters:
546: . svd - singular value solver context obtained from SVDCreate()
548: Output Parameter:
549: . ip - inner product context
551: Level: advanced
553: .seealso: SVDSetIP()
554: @*/
555: PetscErrorCode SVDGetIP(SVD svd,IP *ip)
556: {
560: *ip = svd->ip;
561: return(0);
562: }