1: /*
2: BV (basis vectors) interface routines, callable by users.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, 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_MultVec = 0,BV_MultInPlace = 0,BV_Dot = 0,BV_DotVec = 0,BV_Orthogonalize = 0,BV_OrthogonalizeVec = 0,BV_Scale = 0,BV_Norm = 0,BV_NormVec = 0,BV_SetRandom = 0,BV_MatMult = 0,BV_MatMultVec = 0,BV_MatProject = 0;
28: static PetscBool BVPackageInitialized = PETSC_FALSE;
30: const char *BVOrthogTypes[] = {"CGS","MGS","BVOrthogType","BV_ORTHOG_",0};
31: const char *BVOrthogRefineTypes[] = {"IFNEEDED","NEVER","ALWAYS","BVOrthogRefineType","BV_ORTHOG_REFINE_",0};
32: const char *BVOrthogBlockTypes[] = {"GS","CHOL","BVOrthogBlockType","BV_ORTHOG_BLOCK_",0};
33: const char *BVMatMultTypes[] = {"VECS","MAT","MAT_SAVE","BVMatMultType","BV_MATMULT_",0};
37: /*@C
38: BVFinalizePackage - This function destroys everything in the Slepc interface
39: to the BV package. It is called from SlepcFinalize().
41: Level: developer
43: .seealso: SlepcFinalize()
44: @*/
45: PetscErrorCode BVFinalizePackage(void) 46: {
50: PetscFunctionListDestroy(&BVList);
51: BVPackageInitialized = PETSC_FALSE;
52: BVRegisterAllCalled = PETSC_FALSE;
53: return(0);
54: }
58: /*@C
59: BVInitializePackage - This function initializes everything in the BV package.
60: It is called from PetscDLLibraryRegister() when using dynamic libraries, and
61: on the first call to BVCreate() when using static libraries.
63: Level: developer
65: .seealso: SlepcInitialize()
66: @*/
67: PetscErrorCode BVInitializePackage(void) 68: {
69: char logList[256];
70: char *className;
71: PetscBool opt;
75: if (BVPackageInitialized) return(0);
76: BVPackageInitialized = PETSC_TRUE;
77: /* Register Classes */
78: PetscClassIdRegister("Basis Vectors",&BV_CLASSID);
79: /* Register Constructors */
80: BVRegisterAll();
81: /* Register Events */
82: PetscLogEventRegister("BVCreate",BV_CLASSID,&BV_Create);
83: PetscLogEventRegister("BVCopy",BV_CLASSID,&BV_Copy);
84: PetscLogEventRegister("BVMult",BV_CLASSID,&BV_Mult);
85: PetscLogEventRegister("BVMultVec",BV_CLASSID,&BV_MultVec);
86: PetscLogEventRegister("BVMultInPlace",BV_CLASSID,&BV_MultInPlace);
87: PetscLogEventRegister("BVDot",BV_CLASSID,&BV_Dot);
88: PetscLogEventRegister("BVDotVec",BV_CLASSID,&BV_DotVec);
89: PetscLogEventRegister("BVOrthogonalize",BV_CLASSID,&BV_Orthogonalize);
90: PetscLogEventRegister("BVOrthogonalizeV",BV_CLASSID,&BV_OrthogonalizeVec);
91: PetscLogEventRegister("BVScale",BV_CLASSID,&BV_Scale);
92: PetscLogEventRegister("BVNorm",BV_CLASSID,&BV_Norm);
93: PetscLogEventRegister("BVNormVec",BV_CLASSID,&BV_NormVec);
94: PetscLogEventRegister("BVSetRandom",BV_CLASSID,&BV_SetRandom);
95: PetscLogEventRegister("BVMatMult",BV_CLASSID,&BV_MatMult);
96: PetscLogEventRegister("BVMatMultVec",BV_CLASSID,&BV_MatMultVec);
97: PetscLogEventRegister("BVMatProject",BV_CLASSID,&BV_MatProject);
98: /* Process info exclusions */
99: PetscOptionsGetString(NULL,NULL,"-info_exclude",logList,256,&opt);
100: if (opt) {
101: PetscStrstr(logList,"bv",&className);
102: if (className) {
103: PetscInfoDeactivateClass(BV_CLASSID);
104: }
105: }
106: /* Process summary exclusions */
107: PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,256,&opt);
108: if (opt) {
109: PetscStrstr(logList,"bv",&className);
110: if (className) {
111: PetscLogEventDeactivateClass(BV_CLASSID);
112: }
113: }
114: PetscRegisterFinalize(BVFinalizePackage);
115: return(0);
116: }
120: /*@
121: BVDestroy - Destroys BV context that was created with BVCreate().
123: Collective on BV125: Input Parameter:
126: . bv - the basis vectors context
128: Level: beginner
130: .seealso: BVCreate()
131: @*/
132: PetscErrorCode BVDestroy(BV *bv)133: {
137: if (!*bv) return(0);
139: if (--((PetscObject)(*bv))->refct > 0) { *bv = 0; return(0); }
140: if ((*bv)->ops->destroy) { (*(*bv)->ops->destroy)(*bv); }
141: VecDestroy(&(*bv)->t);
142: MatDestroy(&(*bv)->matrix);
143: VecDestroy(&(*bv)->Bx);
144: BVDestroy(&(*bv)->cached);
145: PetscFree((*bv)->work);
146: PetscFree2((*bv)->h,(*bv)->c);
147: PetscFree((*bv)->omega);
148: MatDestroy(&(*bv)->B);
149: MatDestroy(&(*bv)->C);
150: PetscRandomDestroy(&(*bv)->rand);
151: PetscHeaderDestroy(bv);
152: return(0);
153: }
157: /*@
158: BVCreate - Creates a basis vectors context.
160: Collective on MPI_Comm
162: Input Parameter:
163: . comm - MPI communicator
165: Output Parameter:
166: . bv - location to put the basis vectors context
168: Level: beginner
170: .seealso: BVSetUp(), BVDestroy(), BV171: @*/
172: PetscErrorCode BVCreate(MPI_Comm comm,BV *newbv)173: {
175: BV bv;
179: *newbv = 0;
180: BVInitializePackage();
181: SlepcHeaderCreate(bv,BV_CLASSID,"BV","Basis Vectors","BV",comm,BVDestroy,BVView);
183: bv->t = NULL;
184: bv->n = -1;
185: bv->N = -1;
186: bv->m = 0;
187: bv->l = 0;
188: bv->k = 0;
189: bv->nc = 0;
190: bv->orthog_type = BV_ORTHOG_CGS;
191: bv->orthog_ref = BV_ORTHOG_REFINE_IFNEEDED;
192: bv->orthog_eta = 0.7071;
193: bv->orthog_block = BV_ORTHOG_BLOCK_GS;
194: bv->matrix = NULL;
195: bv->indef = PETSC_FALSE;
196: bv->vmm = BV_MATMULT_MAT;
198: bv->Bx = NULL;
199: bv->xid = 0;
200: bv->xstate = 0;
201: bv->cv[0] = NULL;
202: bv->cv[1] = NULL;
203: bv->ci[0] = -1;
204: bv->ci[1] = -1;
205: bv->st[0] = -1;
206: bv->st[1] = -1;
207: bv->id[0] = 0;
208: bv->id[1] = 0;
209: bv->h = NULL;
210: bv->c = NULL;
211: bv->omega = NULL;
212: bv->B = NULL;
213: bv->C = NULL;
214: bv->Aid = 0;
215: bv->defersfo = PETSC_FALSE;
216: bv->cached = NULL;
217: bv->bvstate = 0;
218: bv->rand = NULL;
219: bv->rrandom = PETSC_FALSE;
220: bv->work = NULL;
221: bv->lwork = 0;
222: bv->data = NULL;
224: *newbv = bv;
225: return(0);
226: }
230: /*@
231: BVInsertVec - Insert a vector into the specified column.
233: Collective on BV235: Input Parameters:
236: + V - basis vectors
237: . j - the column of V to be overwritten
238: - w - the vector to be copied
240: Level: intermediate
242: .seealso: BVInsertVecs()
243: @*/
244: PetscErrorCode BVInsertVec(BV V,PetscInt j,Vec w)245: {
247: PetscInt n,N;
248: Vec v;
255: BVCheckSizes(V,1);
258: VecGetSize(w,&N);
259: VecGetLocalSize(w,&n);
260: 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);
261: 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);
263: BVGetColumn(V,j,&v);
264: VecCopy(w,v);
265: BVRestoreColumn(V,j,&v);
266: PetscObjectStateIncrease((PetscObject)V);
267: return(0);
268: }
272: /*@
273: BVInsertVecs - Insert a set of vectors into the specified columns.
275: Collective on BV277: Input Parameters:
278: + V - basis vectors
279: . s - first column of V to be overwritten
280: . W - set of vectors to be copied
281: - orth - flag indicating if the vectors must be orthogonalized
283: Input/Output Parameter:
284: . m - number of input vectors, on output the number of linearly independent
285: vectors
287: Notes:
288: Copies the contents of vectors W to V(:,s:s+n). If the orthogonalization
289: flag is set, then the vectors are copied one by one and then orthogonalized
290: against the previous ones. If any of them is linearly dependent then it
291: is discarded and the value of m is decreased.
293: Level: intermediate
295: .seealso: BVInsertVec(), BVOrthogonalizeColumn()
296: @*/
297: PetscErrorCode BVInsertVecs(BV V,PetscInt s,PetscInt *m,Vec *W,PetscBool orth)298: {
300: PetscInt n,N,i,ndep;
301: PetscBool lindep;
302: PetscReal norm;
303: Vec v;
310: if (!*m) return(0);
311: if (*m<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",*m);
316: BVCheckSizes(V,1);
319: VecGetSize(*W,&N);
320: VecGetLocalSize(*W,&n);
321: 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);
322: 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);
323: 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);
325: ndep = 0;
326: for (i=0;i<*m;i++) {
327: BVGetColumn(V,s+i-ndep,&v);
328: VecCopy(W[i],v);
329: BVRestoreColumn(V,s+i-ndep,&v);
330: if (orth) {
331: BVOrthogonalizeColumn(V,s+i-ndep,NULL,&norm,&lindep);
332: if (norm==0.0 || lindep) {
333: PetscInfo1(V,"Removing linearly dependent vector %D\n",i);
334: ndep++;
335: } else {
336: BVScaleColumn(V,s+i-ndep,1.0/norm);
337: }
338: }
339: }
340: *m -= ndep;
341: PetscObjectStateIncrease((PetscObject)V);
342: return(0);
343: }
347: /*@
348: BVInsertConstraints - Insert a set of vectors as constraints.
350: Collective on BV352: Input Parameters:
353: + V - basis vectors
354: - C - set of vectors to be inserted as constraints
356: Input/Output Parameter:
357: . nc - number of input vectors, on output the number of linearly independent
358: vectors
360: Notes:
361: The constraints are relevant only during orthogonalization. Constraint
362: vectors span a subspace that is deflated in every orthogonalization
363: operation, so they are intended for removing those directions from the
364: orthogonal basis computed in regular BV columns.
366: Constraints are not stored in regular BV colums, but in a special part of
367: the storage. They can be accessed with negative indices in BVGetColumn().
369: This operation is DESTRUCTIVE, meaning that all data contained in the
370: columns of V is lost. This is typically invoked just after creating the BV.
371: Once a set of constraints has been set, it is not allowed to call this
372: function again.
374: The vectors are copied one by one and then orthogonalized against the
375: previous ones. If any of them is linearly dependent then it is discarded
376: and the value of nc is decreased. The behaviour is similar to BVInsertVecs().
378: Level: advanced
380: .seealso: BVInsertVecs(), BVOrthogonalizeColumn(), BVGetColumn(), BVGetNumConstraints()
381: @*/
382: PetscErrorCode BVInsertConstraints(BV V,PetscInt *nc,Vec *C)383: {
385: PetscInt msave;
391: if (!*nc) return(0);
392: if (*nc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of constraints (given %D) cannot be negative",*nc);
396: BVCheckSizes(V,1);
398: if (V->nc) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_WRONGSTATE,"Constraints already present in this BV object");
399: if (V->ci[0]!=-1 || V->ci[1]!=-1) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Cannot call BVInsertConstraints after BVGetColumn");
401: msave = V->m;
402: BVResize(V,*nc+V->m,PETSC_FALSE);
403: BVInsertVecs(V,0,nc,C,PETSC_TRUE);
404: V->nc = *nc;
405: V->m = msave;
406: V->ci[0] = -V->nc-1;
407: V->ci[1] = -V->nc-1;
408: PetscObjectStateIncrease((PetscObject)V);
409: return(0);
410: }
414: /*@C
415: BVSetOptionsPrefix - Sets the prefix used for searching for all
416: BV options in the database.
418: Logically Collective on BV420: Input Parameters:
421: + bv - the basis vectors context
422: - prefix - the prefix string to prepend to all BV option requests
424: Notes:
425: A hyphen (-) must NOT be given at the beginning of the prefix name.
426: The first character of all runtime options is AUTOMATICALLY the
427: hyphen.
429: Level: advanced
431: .seealso: BVAppendOptionsPrefix(), BVGetOptionsPrefix()
432: @*/
433: PetscErrorCode BVSetOptionsPrefix(BV bv,const char *prefix)434: {
439: PetscObjectSetOptionsPrefix((PetscObject)bv,prefix);
440: return(0);
441: }
445: /*@C
446: BVAppendOptionsPrefix - Appends to the prefix used for searching for all
447: BV options in the database.
449: Logically Collective on BV451: Input Parameters:
452: + bv - the basis vectors context
453: - prefix - the prefix string to prepend to all BV option requests
455: Notes:
456: A hyphen (-) must NOT be given at the beginning of the prefix name.
457: The first character of all runtime options is AUTOMATICALLY the
458: hyphen.
460: Level: advanced
462: .seealso: BVSetOptionsPrefix(), BVGetOptionsPrefix()
463: @*/
464: PetscErrorCode BVAppendOptionsPrefix(BV bv,const char *prefix)465: {
470: PetscObjectAppendOptionsPrefix((PetscObject)bv,prefix);
471: return(0);
472: }
476: /*@C
477: BVGetOptionsPrefix - Gets the prefix used for searching for all
478: BV options in the database.
480: Not Collective
482: Input Parameters:
483: . bv - the basis vectors context
485: Output Parameters:
486: . prefix - pointer to the prefix string used, is returned
488: Note:
489: On the Fortran side, the user should pass in a string 'prefix' of
490: sufficient length to hold the prefix.
492: Level: advanced
494: .seealso: BVSetOptionsPrefix(), BVAppendOptionsPrefix()
495: @*/
496: PetscErrorCode BVGetOptionsPrefix(BV bv,const char *prefix[])497: {
503: PetscObjectGetOptionsPrefix((PetscObject)bv,prefix);
504: return(0);
505: }
509: static PetscErrorCode BVView_Default(BV bv,PetscViewer viewer)510: {
511: PetscErrorCode ierr;
512: PetscInt j;
513: Vec v;
514: PetscViewerFormat format;
515: PetscBool isascii,ismatlab=PETSC_FALSE;
516: const char *bvname,*name;
519: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
520: if (isascii) {
521: PetscViewerGetFormat(viewer,&format);
522: if (format == PETSC_VIEWER_ASCII_MATLAB) ismatlab = PETSC_TRUE;
523: }
524: if (ismatlab) {
525: PetscObjectGetName((PetscObject)bv,&bvname);
526: PetscViewerASCIIPrintf(viewer,"%s=[];\n",bvname);
527: }
528: for (j=-bv->nc;j<bv->m;j++) {
529: BVGetColumn(bv,j,&v);
530: VecView(v,viewer);
531: if (ismatlab) {
532: PetscObjectGetName((PetscObject)v,&name);
533: PetscViewerASCIIPrintf(viewer,"%s=[%s,%s];clear %s\n",bvname,bvname,name,name);
534: }
535: BVRestoreColumn(bv,j,&v);
536: }
537: return(0);
538: }
542: /*@C
543: BVView - Prints the BV data structure.
545: Collective on BV547: Input Parameters:
548: + bv - the BV context
549: - viewer - optional visualization context
551: Note:
552: The available visualization contexts include
553: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
554: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
555: output where only the first processor opens
556: the file. All other processors send their
557: data to the first processor to print.
559: The user can open an alternative visualization contexts with
560: PetscViewerASCIIOpen() (output to a specified file).
562: Level: beginner
564: .seealso: PetscViewerASCIIOpen()
565: @*/
566: PetscErrorCode BVView(BV bv,PetscViewer viewer)567: {
568: PetscErrorCode ierr;
569: PetscBool isascii;
570: PetscViewerFormat format;
571: const char *orthname[2] = {"classical","modified"};
572: const char *refname[3] = {"if needed","never","always"};
573: const char *borthname[2] = {"Gram-Schmidt","Cholesky"};
577: if (!viewer) {
578: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)bv),&viewer);
579: }
582: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
583: if (isascii) {
584: PetscObjectPrintClassNamePrefixType((PetscObject)bv,viewer);
585: PetscViewerGetFormat(viewer,&format);
586: PetscViewerASCIIPushTab(viewer);
587: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
588: PetscViewerASCIIPrintf(viewer,"%D columns of global length %D\n",bv->m,bv->N);
589: if (bv->nc>0) {
590: PetscViewerASCIIPrintf(viewer,"number of constraints: %D\n",bv->nc);
591: }
592: PetscViewerASCIIPrintf(viewer,"vector orthogonalization method: %s Gram-Schmidt\n",orthname[bv->orthog_type]);
593: switch (bv->orthog_ref) {
594: case BV_ORTHOG_REFINE_IFNEEDED:
595: PetscViewerASCIIPrintf(viewer,"orthogonalization refinement: %s (eta: %g)\n",refname[bv->orthog_ref],(double)bv->orthog_eta);
596: break;
597: case BV_ORTHOG_REFINE_NEVER:
598: case BV_ORTHOG_REFINE_ALWAYS:
599: PetscViewerASCIIPrintf(viewer,"orthogonalization refinement: %s\n",refname[bv->orthog_ref]);
600: break;
601: }
602: PetscViewerASCIIPrintf(viewer,"block orthogonalization method: %s\n",borthname[bv->orthog_block]);
603: if (bv->matrix) {
604: if (bv->indef) {
605: PetscViewerASCIIPrintf(viewer,"indefinite inner product\n");
606: } else {
607: PetscViewerASCIIPrintf(viewer,"non-standard inner product\n");
608: }
609: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
610: MatView(bv->matrix,viewer);
611: PetscViewerPopFormat(viewer);
612: }
613: switch (bv->vmm) {
614: case BV_MATMULT_VECS:
615: PetscViewerASCIIPrintf(viewer,"doing matmult as matrix-vector products\n");
616: break;
617: case BV_MATMULT_MAT:
618: PetscViewerASCIIPrintf(viewer,"doing matmult as a single matrix-matrix product\n");
619: break;
620: case BV_MATMULT_MAT_SAVE:
621: PetscViewerASCIIPrintf(viewer,"doing matmult as a single matrix-matrix product, saving aux matrices\n");
622: break;
623: }
624: if (bv->rrandom) {
625: PetscViewerASCIIPrintf(viewer,"generating random vectors independent of the number of processes\n");
626: }
627: } else {
628: if (bv->ops->view) { (*bv->ops->view)(bv,viewer); }
629: else { BVView_Default(bv,viewer); }
630: }
631: PetscViewerASCIIPopTab(viewer);
632: } else {
633: (*bv->ops->view)(bv,viewer);
634: }
635: return(0);
636: }
640: /*@C
641: BVRegister - Adds a new storage format to de BV package.
643: Not collective
645: Input Parameters:
646: + name - name of a new user-defined BV647: - function - routine to create context
649: Notes:
650: BVRegister() may be called multiple times to add several user-defined
651: basis vectors.
653: Level: advanced
655: .seealso: BVRegisterAll()
656: @*/
657: PetscErrorCode BVRegister(const char *name,PetscErrorCode (*function)(BV))658: {
662: PetscFunctionListAdd(&BVList,name,function);
663: return(0);
664: }
668: PetscErrorCode BVAllocateWork_Private(BV bv,PetscInt s)669: {
673: if (s>bv->lwork) {
674: PetscFree(bv->work);
675: PetscMalloc1(s,&bv->work);
676: PetscLogObjectMemory((PetscObject)bv,(s-bv->lwork)*sizeof(PetscScalar));
677: bv->lwork = s;
678: }
679: return(0);
680: }