1: /*
2: BV (basis vectors) interface routines, callable by users.
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/bvimpl.h> /*I "slepcbv.h" I*/
26: PetscClassId BV_CLASSID = 0;
27: PetscLogEvent BV_Create = 0,BV_Copy = 0,BV_Mult = 0,BV_Dot = 0,BV_Orthogonalize = 0,BV_Scale = 0,BV_Norm = 0,BV_SetRandom = 0,BV_MatMult = 0,BV_MatProject = 0,BV_AXPY = 0;
28: static PetscBool BVPackageInitialized = PETSC_FALSE;
32: /*@C
33: BVFinalizePackage - This function destroys everything in the Slepc interface
34: to the BV package. It is called from SlepcFinalize().
36: Level: developer
38: .seealso: SlepcFinalize()
39: @*/
40: PetscErrorCode BVFinalizePackage(void) 41: {
45: PetscFunctionListDestroy(&BVList);
46: BVPackageInitialized = PETSC_FALSE;
47: BVRegisterAllCalled = PETSC_FALSE;
48: return(0);
49: }
53: /*@C
54: BVInitializePackage - This function initializes everything in the BV package.
55: It is called from PetscDLLibraryRegister() when using dynamic libraries, and
56: on the first call to BVCreate() when using static libraries.
58: Level: developer
60: .seealso: SlepcInitialize()
61: @*/
62: PetscErrorCode BVInitializePackage(void) 63: {
64: char logList[256];
65: char *className;
66: PetscBool opt;
70: if (BVPackageInitialized) return(0);
71: BVPackageInitialized = PETSC_TRUE;
72: /* Register Classes */
73: PetscClassIdRegister("Basis Vectors",&BV_CLASSID);
74: /* Register Constructors */
75: BVRegisterAll();
76: /* Register Events */
77: PetscLogEventRegister("BVCreate",BV_CLASSID,&BV_Create);
78: PetscLogEventRegister("BVCopy",BV_CLASSID,&BV_Copy);
79: PetscLogEventRegister("BVMult",BV_CLASSID,&BV_Mult);
80: PetscLogEventRegister("BVDot",BV_CLASSID,&BV_Dot);
81: PetscLogEventRegister("BVOrthogonalize",BV_CLASSID,&BV_Orthogonalize);
82: PetscLogEventRegister("BVScale",BV_CLASSID,&BV_Scale);
83: PetscLogEventRegister("BVNorm",BV_CLASSID,&BV_Norm);
84: PetscLogEventRegister("BVSetRandom",BV_CLASSID,&BV_SetRandom);
85: PetscLogEventRegister("BVMatMult",BV_CLASSID,&BV_MatMult);
86: PetscLogEventRegister("BVMatProject",BV_CLASSID,&BV_MatProject);
87: PetscLogEventRegister("BVAXPY",BV_CLASSID,&BV_AXPY);
88: /* Process info exclusions */
89: PetscOptionsGetString(NULL,"-info_exclude",logList,256,&opt);
90: if (opt) {
91: PetscStrstr(logList,"bv",&className);
92: if (className) {
93: PetscInfoDeactivateClass(BV_CLASSID);
94: }
95: }
96: /* Process summary exclusions */
97: PetscOptionsGetString(NULL,"-log_summary_exclude",logList,256,&opt);
98: if (opt) {
99: PetscStrstr(logList,"bv",&className);
100: if (className) {
101: PetscLogEventDeactivateClass(BV_CLASSID);
102: }
103: }
104: PetscRegisterFinalize(BVFinalizePackage);
105: return(0);
106: }
110: /*@C
111: BVDestroy - Destroys BV context that was created with BVCreate().
113: Collective on BV115: Input Parameter:
116: . bv - the basis vectors context
118: Level: beginner
120: .seealso: BVCreate()
121: @*/
122: PetscErrorCode BVDestroy(BV *bv)123: {
127: if (!*bv) return(0);
129: if (--((PetscObject)(*bv))->refct > 0) { *bv = 0; return(0); }
130: if ((*bv)->ops->destroy) { (*(*bv)->ops->destroy)(*bv); }
131: VecDestroy(&(*bv)->t);
132: MatDestroy(&(*bv)->matrix);
133: VecDestroy(&(*bv)->Bx);
134: PetscFree((*bv)->work);
135: PetscFree2((*bv)->h,(*bv)->c);
136: PetscFree((*bv)->omega);
137: PetscHeaderDestroy(bv);
138: return(0);
139: }
143: /*@C
144: BVCreate - Creates a basis vectors context.
146: Collective on MPI_Comm
148: Input Parameter:
149: . comm - MPI communicator
151: Output Parameter:
152: . bv - location to put the basis vectors context
154: Level: beginner
156: .seealso: BVSetUp(), BVDestroy(), BV157: @*/
158: PetscErrorCode BVCreate(MPI_Comm comm,BV *newbv)159: {
161: BV bv;
165: *newbv = 0;
166: BVInitializePackage();
167: SlepcHeaderCreate(bv,_p_BV,struct _BVOps,BV_CLASSID,"BV","Basis Vectors","BV",comm,BVDestroy,BVView);
169: bv->t = NULL;
170: bv->n = -1;
171: bv->N = -1;
172: bv->m = 0;
173: bv->l = 0;
174: bv->k = 0;
175: bv->nc = 0;
176: bv->orthog_type = BV_ORTHOG_CGS;
177: bv->orthog_ref = BV_ORTHOG_REFINE_IFNEEDED;
178: bv->orthog_eta = 0.7071;
179: bv->matrix = NULL;
180: bv->indef = PETSC_FALSE;
182: bv->Bx = NULL;
183: bv->xid = 0;
184: bv->xstate = 0;
185: bv->cv[0] = NULL;
186: bv->cv[1] = NULL;
187: bv->ci[0] = -1;
188: bv->ci[1] = -1;
189: bv->st[0] = -1;
190: bv->st[1] = -1;
191: bv->id[0] = 0;
192: bv->id[1] = 0;
193: bv->h = NULL;
194: bv->c = NULL;
195: bv->omega = NULL;
196: bv->work = NULL;
197: bv->lwork = 0;
198: bv->data = NULL;
200: *newbv = bv;
201: return(0);
202: }
206: /*@
207: BVInsertVec - Insert a vector into the specified column.
209: Collective on BV211: Input Parameters:
212: + V - basis vectors
213: . j - the column of V to be overwritten
214: - w - the vector to be copied
216: Level: intermediate
218: .seealso: BVInsertVecs()
219: @*/
220: PetscErrorCode BVInsertVec(BV V,PetscInt j,Vec w)221: {
223: PetscInt n,N;
224: Vec v;
231: BVCheckSizes(V,1);
234: VecGetSize(w,&N);
235: VecGetLocalSize(w,&n);
236: if (N!=V->N || n!=V->n) SETERRQ4(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Vec sizes (global %D, local %D) do not match BV sizes (global %D, local %D)",N,n,V->N,V->n);
237: if (j<-V->nc || j>=V->m) SETERRQ3(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument j has wrong value %D, should be between %D and %D",j,-V->nc,V->m-1);
239: BVGetColumn(V,j,&v);
240: VecCopy(w,v);
241: BVRestoreColumn(V,j,&v);
242: PetscObjectStateIncrease((PetscObject)V);
243: return(0);
244: }
248: /*@
249: BVInsertVecs - Insert a set of vectors into the specified columns.
251: Collective on BV253: Input Parameters:
254: + V - basis vectors
255: . s - first column of V to be overwritten
256: . W - set of vectors to be copied
257: - orth - flag indicating if the vectors must be orthogonalized
259: Input/Output Parameter:
260: . m - number of input vectors, on output the number of linearly independent
261: vectors
263: Notes:
264: Copies the contents of vectors W to V(:,s:s+n). If the orthogonalization
265: flag is set, then the vectors are copied one by one and then orthogonalized
266: against the previous ones. If any of them is linearly dependent then it
267: is discarded and the value of m is decreased.
269: Level: intermediate
271: .seealso: BVInsertVec(), BVOrthogonalizeColumn()
272: @*/
273: PetscErrorCode BVInsertVecs(BV V,PetscInt s,PetscInt *m,Vec *W,PetscBool orth)274: {
276: PetscInt n,N,i,ndep;
277: PetscBool lindep;
278: PetscReal norm;
279: Vec v;
286: if (!*m) return(0);
287: if (*m<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",*m);
292: BVCheckSizes(V,1);
295: VecGetSize(*W,&N);
296: VecGetLocalSize(*W,&n);
297: if (N!=V->N || n!=V->n) SETERRQ4(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Vec sizes (global %D, local %D) do not match BV sizes (global %D, local %D)",N,n,V->N,V->n);
298: if (s<0 || s>=V->m) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument s has wrong value %D, should be between 0 and %D",s,V->m-1);
299: if (s+(*m)>V->m) SETERRQ1(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Too many vectors provided, there is only room for %D",V->m);
301: ndep = 0;
302: for (i=0;i<*m;i++) {
303: BVGetColumn(V,s+i-ndep,&v);
304: VecCopy(W[i],v);
305: BVRestoreColumn(V,s+i-ndep,&v);
306: if (orth) {
307: BVOrthogonalizeColumn(V,s+i-ndep,NULL,&norm,&lindep);
308: if (norm==0.0 || lindep) {
309: PetscInfo1(V,"Removing linearly dependent vector %D\n",i);
310: ndep++;
311: } else {
312: BVScaleColumn(V,s+i-ndep,1.0/norm);
313: }
314: }
315: }
316: *m -= ndep;
317: PetscObjectStateIncrease((PetscObject)V);
318: return(0);
319: }
323: /*@
324: BVInsertConstraints - Insert a set of vectors as constraints.
326: Collective on BV328: Input Parameters:
329: + V - basis vectors
330: - C - set of vectors to be inserted as constraints
332: Input/Output Parameter:
333: . nc - number of input vectors, on output the number of linearly independent
334: vectors
336: Notes:
337: The constraints are relevant only during orthogonalization. Constraint
338: vectors span a subspace that is deflated in every orthogonalization
339: operation, so they are intended for removing those directions from the
340: orthogonal basis computed in regular BV columns.
342: Constraints are not stored in regular BV colums, but in a special part of
343: the storage. They can be accessed with negative indices in BVGetColumn().
345: This operation is DESTRUCTIVE, meaning that all data contained in the
346: columns of V is lost. This is typically invoked just after creating the BV.
347: Once a set of constraints has been set, it is not allowed to call this
348: function again.
350: The vectors are copied one by one and then orthogonalized against the
351: previous ones. If any of them is linearly dependent then it is discarded
352: and the value of nc is decreased. The behaviour is similar to BVInsertVecs().
354: Level: advanced
356: .seealso: BVInsertVecs(), BVOrthogonalizeColumn(), BVGetColumn(), BVGetNumConstraints()
357: @*/
358: PetscErrorCode BVInsertConstraints(BV V,PetscInt *nc,Vec *C)359: {
361: PetscInt msave;
367: if (!*nc) return(0);
368: if (*nc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of constraints (given %D) cannot be negative",*nc);
372: BVCheckSizes(V,1);
374: if (V->nc) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_WRONGSTATE,"Constraints already present in this BV object");
375: if (V->ci[0]!=-1 || V->ci[1]!=-1) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Cannot call BVInsertConstraints after BVGetColumn");
377: msave = V->m;
378: BVResize(V,*nc+V->m,PETSC_FALSE);
379: BVInsertVecs(V,0,nc,C,PETSC_TRUE);
380: V->nc = *nc;
381: V->m = msave;
382: V->ci[0] = -V->nc-1;
383: V->ci[1] = -V->nc-1;
384: PetscObjectStateIncrease((PetscObject)V);
385: return(0);
386: }
390: /*@C
391: BVSetOptionsPrefix - Sets the prefix used for searching for all
392: BV options in the database.
394: Logically Collective on BV396: Input Parameters:
397: + bv - the basis vectors context
398: - prefix - the prefix string to prepend to all BV option requests
400: Notes:
401: A hyphen (-) must NOT be given at the beginning of the prefix name.
402: The first character of all runtime options is AUTOMATICALLY the
403: hyphen.
405: Level: advanced
407: .seealso: BVAppendOptionsPrefix(), BVGetOptionsPrefix()
408: @*/
409: PetscErrorCode BVSetOptionsPrefix(BV bv,const char *prefix)410: {
415: PetscObjectSetOptionsPrefix((PetscObject)bv,prefix);
416: return(0);
417: }
421: /*@C
422: BVAppendOptionsPrefix - Appends to the prefix used for searching for all
423: BV options in the database.
425: Logically Collective on BV427: Input Parameters:
428: + bv - the basis vectors context
429: - prefix - the prefix string to prepend to all BV option requests
431: Notes:
432: A hyphen (-) must NOT be given at the beginning of the prefix name.
433: The first character of all runtime options is AUTOMATICALLY the
434: hyphen.
436: Level: advanced
438: .seealso: BVSetOptionsPrefix(), BVGetOptionsPrefix()
439: @*/
440: PetscErrorCode BVAppendOptionsPrefix(BV bv,const char *prefix)441: {
446: PetscObjectAppendOptionsPrefix((PetscObject)bv,prefix);
447: return(0);
448: }
452: /*@C
453: BVGetOptionsPrefix - Gets the prefix used for searching for all
454: BV options in the database.
456: Not Collective
458: Input Parameters:
459: . bv - the basis vectors context
461: Output Parameters:
462: . prefix - pointer to the prefix string used, is returned
464: Notes: On the Fortran side, the user should pass in a string 'prefix' of
465: sufficient length to hold the prefix.
467: Level: advanced
469: .seealso: BVSetOptionsPrefix(), BVAppendOptionsPrefix()
470: @*/
471: PetscErrorCode BVGetOptionsPrefix(BV bv,const char *prefix[])472: {
478: PetscObjectGetOptionsPrefix((PetscObject)bv,prefix);
479: return(0);
480: }
484: static PetscErrorCode BVView_Default(BV bv,PetscViewer viewer)485: {
486: PetscErrorCode ierr;
487: PetscInt j;
488: Vec v;
489: PetscViewerFormat format;
490: PetscBool isascii,ismatlab=PETSC_FALSE;
493: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
494: if (isascii) {
495: PetscViewerGetFormat(viewer,&format);
496: if (format == PETSC_VIEWER_ASCII_MATLAB) ismatlab = PETSC_TRUE;
497: }
498: if (ismatlab) {
499: PetscViewerASCIIPrintf(viewer,"%s=[];\n",((PetscObject)bv)->name);
500: }
501: for (j=bv->nc;j<bv->nc+bv->m;j++) {
502: BVGetColumn(bv,j,&v);
503: VecView(v,viewer);
504: if (ismatlab) {
505: PetscViewerASCIIPrintf(viewer,"%s=[%s,%s];clear %s\n",((PetscObject)bv)->name,((PetscObject)bv)->name,((PetscObject)v)->name,((PetscObject)v)->name);
506: }
507: BVRestoreColumn(bv,j,&v);
508: }
509: return(0);
510: }
514: /*@C
515: BVView - Prints the BV data structure.
517: Collective on BV519: Input Parameters:
520: + bv - the BV context
521: - viewer - optional visualization context
523: Note:
524: The available visualization contexts include
525: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
526: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
527: output where only the first processor opens
528: the file. All other processors send their
529: data to the first processor to print.
531: The user can open an alternative visualization contexts with
532: PetscViewerASCIIOpen() (output to a specified file).
534: Level: beginner
536: .seealso: PetscViewerASCIIOpen()
537: @*/
538: PetscErrorCode BVView(BV bv,PetscViewer viewer)539: {
540: PetscErrorCode ierr;
541: PetscBool isascii;
542: PetscViewerFormat format;
543: const char *orthname[2] = {"classical","modified"};
544: const char *refname[3] = {"if needed","never","always"};
549: if (!viewer) {
550: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)bv),&viewer);
551: }
554: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
555: if (isascii) {
556: PetscObjectPrintClassNamePrefixType((PetscObject)bv,viewer);
557: PetscViewerGetFormat(viewer,&format);
558: PetscViewerASCIIPushTab(viewer);
559: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
560: PetscViewerASCIIPrintf(viewer,"%D columns of global length %D\n",bv->m,bv->N);
561: if (bv->nc>0) {
562: PetscViewerASCIIPrintf(viewer,"number of constraints: %D\n",bv->nc);
563: }
564: PetscViewerASCIIPrintf(viewer,"orthogonalization method: %s Gram-Schmidt\n",orthname[bv->orthog_type]);
565: switch (bv->orthog_ref) {
566: case BV_ORTHOG_REFINE_IFNEEDED:
567: PetscViewerASCIIPrintf(viewer,"orthogonalization refinement: %s (eta: %g)\n",refname[bv->orthog_ref],(double)bv->orthog_eta);
568: break;
569: case BV_ORTHOG_REFINE_NEVER:
570: case BV_ORTHOG_REFINE_ALWAYS:
571: PetscViewerASCIIPrintf(viewer,"orthogonalization refinement: %s\n",refname[bv->orthog_ref]);
572: break;
573: }
574: if (bv->matrix) {
575: if (bv->indef) {
576: PetscViewerASCIIPrintf(viewer,"indefinite inner product\n");
577: } else {
578: PetscViewerASCIIPrintf(viewer,"non-standard inner product\n");
579: }
580: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
581: MatView(bv->matrix,viewer);
582: PetscViewerPopFormat(viewer);
583: }
584: } else {
585: if (bv->ops->view) { (*bv->ops->view)(bv,viewer); }
586: else { BVView_Default(bv,viewer); }
587: }
588: PetscViewerASCIIPopTab(viewer);
589: } else {
590: (*bv->ops->view)(bv,viewer);
591: }
592: return(0);
593: }
597: /*@C
598: BVRegister - Adds a new storage format to de BV package.
600: Not collective
602: Input Parameters:
603: + name - name of a new user-defined BV604: - function - routine to create context
606: Notes:
607: BVRegister() may be called multiple times to add several user-defined
608: basis vectors.
610: Level: advanced
612: .seealso: BVRegisterAll()
613: @*/
614: PetscErrorCode BVRegister(const char *name,PetscErrorCode (*function)(BV))615: {
619: PetscFunctionListAdd(&BVList,name,function);
620: return(0);
621: }
625: PetscErrorCode BVAllocateWork_Private(BV bv,PetscInt s)626: {
630: if (s>bv->lwork) {
631: PetscFree(bv->work);
632: PetscMalloc1(s,&bv->work);
633: PetscLogObjectMemory((PetscObject)bv,(s-bv->lwork)*sizeof(PetscScalar));
634: bv->lwork = s;
635: }
636: return(0);
637: }