1: /*
2: BV (basis vectors) interface routines, callable by users.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, 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: /*@
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: BVDestroy(&(*bv)->cached);
135: PetscFree((*bv)->work);
136: PetscFree2((*bv)->h,(*bv)->c);
137: PetscFree((*bv)->omega);
138: MatDestroy(&(*bv)->B);
139: MatDestroy(&(*bv)->C);
140: PetscHeaderDestroy(bv);
141: return(0);
142: }
146: /*@
147: BVCreate - Creates a basis vectors context.
149: Collective on MPI_Comm
151: Input Parameter:
152: . comm - MPI communicator
154: Output Parameter:
155: . bv - location to put the basis vectors context
157: Level: beginner
159: .seealso: BVSetUp(), BVDestroy(), BV160: @*/
161: PetscErrorCode BVCreate(MPI_Comm comm,BV *newbv)162: {
164: BV bv;
168: *newbv = 0;
169: BVInitializePackage();
170: SlepcHeaderCreate(bv,BV_CLASSID,"BV","Basis Vectors","BV",comm,BVDestroy,BVView);
172: bv->t = NULL;
173: bv->n = -1;
174: bv->N = -1;
175: bv->m = 0;
176: bv->l = 0;
177: bv->k = 0;
178: bv->nc = 0;
179: bv->orthog_type = BV_ORTHOG_CGS;
180: bv->orthog_ref = BV_ORTHOG_REFINE_IFNEEDED;
181: bv->orthog_eta = 0.7071;
182: bv->orthog_block = BV_ORTHOG_BLOCK_GS;
183: bv->matrix = NULL;
184: bv->indef = PETSC_FALSE;
185: bv->vmm = BV_MATMULT_MAT;
187: bv->Bx = NULL;
188: bv->xid = 0;
189: bv->xstate = 0;
190: bv->cv[0] = NULL;
191: bv->cv[1] = NULL;
192: bv->ci[0] = -1;
193: bv->ci[1] = -1;
194: bv->st[0] = -1;
195: bv->st[1] = -1;
196: bv->id[0] = 0;
197: bv->id[1] = 0;
198: bv->h = NULL;
199: bv->c = NULL;
200: bv->omega = NULL;
201: bv->B = NULL;
202: bv->C = NULL;
203: bv->Aid = 0;
204: bv->defersfo = PETSC_FALSE;
205: bv->cached = NULL;
206: bv->bvstate = 0;
207: bv->work = NULL;
208: bv->lwork = 0;
209: bv->data = NULL;
211: *newbv = bv;
212: return(0);
213: }
217: /*@
218: BVInsertVec - Insert a vector into the specified column.
220: Collective on BV222: Input Parameters:
223: + V - basis vectors
224: . j - the column of V to be overwritten
225: - w - the vector to be copied
227: Level: intermediate
229: .seealso: BVInsertVecs()
230: @*/
231: PetscErrorCode BVInsertVec(BV V,PetscInt j,Vec w)232: {
234: PetscInt n,N;
235: Vec v;
242: BVCheckSizes(V,1);
245: VecGetSize(w,&N);
246: VecGetLocalSize(w,&n);
247: 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);
248: 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);
250: BVGetColumn(V,j,&v);
251: VecCopy(w,v);
252: BVRestoreColumn(V,j,&v);
253: PetscObjectStateIncrease((PetscObject)V);
254: return(0);
255: }
259: /*@
260: BVInsertVecs - Insert a set of vectors into the specified columns.
262: Collective on BV264: Input Parameters:
265: + V - basis vectors
266: . s - first column of V to be overwritten
267: . W - set of vectors to be copied
268: - orth - flag indicating if the vectors must be orthogonalized
270: Input/Output Parameter:
271: . m - number of input vectors, on output the number of linearly independent
272: vectors
274: Notes:
275: Copies the contents of vectors W to V(:,s:s+n). If the orthogonalization
276: flag is set, then the vectors are copied one by one and then orthogonalized
277: against the previous ones. If any of them is linearly dependent then it
278: is discarded and the value of m is decreased.
280: Level: intermediate
282: .seealso: BVInsertVec(), BVOrthogonalizeColumn()
283: @*/
284: PetscErrorCode BVInsertVecs(BV V,PetscInt s,PetscInt *m,Vec *W,PetscBool orth)285: {
287: PetscInt n,N,i,ndep;
288: PetscBool lindep;
289: PetscReal norm;
290: Vec v;
297: if (!*m) return(0);
298: if (*m<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",*m);
303: BVCheckSizes(V,1);
306: VecGetSize(*W,&N);
307: VecGetLocalSize(*W,&n);
308: 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);
309: 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);
310: 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);
312: ndep = 0;
313: for (i=0;i<*m;i++) {
314: BVGetColumn(V,s+i-ndep,&v);
315: VecCopy(W[i],v);
316: BVRestoreColumn(V,s+i-ndep,&v);
317: if (orth) {
318: BVOrthogonalizeColumn(V,s+i-ndep,NULL,&norm,&lindep);
319: if (norm==0.0 || lindep) {
320: PetscInfo1(V,"Removing linearly dependent vector %D\n",i);
321: ndep++;
322: } else {
323: BVScaleColumn(V,s+i-ndep,1.0/norm);
324: }
325: }
326: }
327: *m -= ndep;
328: PetscObjectStateIncrease((PetscObject)V);
329: return(0);
330: }
334: /*@
335: BVInsertConstraints - Insert a set of vectors as constraints.
337: Collective on BV339: Input Parameters:
340: + V - basis vectors
341: - C - set of vectors to be inserted as constraints
343: Input/Output Parameter:
344: . nc - number of input vectors, on output the number of linearly independent
345: vectors
347: Notes:
348: The constraints are relevant only during orthogonalization. Constraint
349: vectors span a subspace that is deflated in every orthogonalization
350: operation, so they are intended for removing those directions from the
351: orthogonal basis computed in regular BV columns.
353: Constraints are not stored in regular BV colums, but in a special part of
354: the storage. They can be accessed with negative indices in BVGetColumn().
356: This operation is DESTRUCTIVE, meaning that all data contained in the
357: columns of V is lost. This is typically invoked just after creating the BV.
358: Once a set of constraints has been set, it is not allowed to call this
359: function again.
361: The vectors are copied one by one and then orthogonalized against the
362: previous ones. If any of them is linearly dependent then it is discarded
363: and the value of nc is decreased. The behaviour is similar to BVInsertVecs().
365: Level: advanced
367: .seealso: BVInsertVecs(), BVOrthogonalizeColumn(), BVGetColumn(), BVGetNumConstraints()
368: @*/
369: PetscErrorCode BVInsertConstraints(BV V,PetscInt *nc,Vec *C)370: {
372: PetscInt msave;
378: if (!*nc) return(0);
379: if (*nc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of constraints (given %D) cannot be negative",*nc);
383: BVCheckSizes(V,1);
385: if (V->nc) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_WRONGSTATE,"Constraints already present in this BV object");
386: if (V->ci[0]!=-1 || V->ci[1]!=-1) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Cannot call BVInsertConstraints after BVGetColumn");
388: msave = V->m;
389: BVResize(V,*nc+V->m,PETSC_FALSE);
390: BVInsertVecs(V,0,nc,C,PETSC_TRUE);
391: V->nc = *nc;
392: V->m = msave;
393: V->ci[0] = -V->nc-1;
394: V->ci[1] = -V->nc-1;
395: PetscObjectStateIncrease((PetscObject)V);
396: return(0);
397: }
401: /*@C
402: BVSetOptionsPrefix - Sets the prefix used for searching for all
403: BV options in the database.
405: Logically Collective on BV407: Input Parameters:
408: + bv - the basis vectors context
409: - prefix - the prefix string to prepend to all BV option requests
411: Notes:
412: A hyphen (-) must NOT be given at the beginning of the prefix name.
413: The first character of all runtime options is AUTOMATICALLY the
414: hyphen.
416: Level: advanced
418: .seealso: BVAppendOptionsPrefix(), BVGetOptionsPrefix()
419: @*/
420: PetscErrorCode BVSetOptionsPrefix(BV bv,const char *prefix)421: {
426: PetscObjectSetOptionsPrefix((PetscObject)bv,prefix);
427: return(0);
428: }
432: /*@C
433: BVAppendOptionsPrefix - Appends to the prefix used for searching for all
434: BV options in the database.
436: Logically Collective on BV438: Input Parameters:
439: + bv - the basis vectors context
440: - prefix - the prefix string to prepend to all BV option requests
442: Notes:
443: A hyphen (-) must NOT be given at the beginning of the prefix name.
444: The first character of all runtime options is AUTOMATICALLY the
445: hyphen.
447: Level: advanced
449: .seealso: BVSetOptionsPrefix(), BVGetOptionsPrefix()
450: @*/
451: PetscErrorCode BVAppendOptionsPrefix(BV bv,const char *prefix)452: {
457: PetscObjectAppendOptionsPrefix((PetscObject)bv,prefix);
458: return(0);
459: }
463: /*@C
464: BVGetOptionsPrefix - Gets the prefix used for searching for all
465: BV options in the database.
467: Not Collective
469: Input Parameters:
470: . bv - the basis vectors context
472: Output Parameters:
473: . prefix - pointer to the prefix string used, is returned
475: Notes: On the Fortran side, the user should pass in a string 'prefix' of
476: sufficient length to hold the prefix.
478: Level: advanced
480: .seealso: BVSetOptionsPrefix(), BVAppendOptionsPrefix()
481: @*/
482: PetscErrorCode BVGetOptionsPrefix(BV bv,const char *prefix[])483: {
489: PetscObjectGetOptionsPrefix((PetscObject)bv,prefix);
490: return(0);
491: }
495: static PetscErrorCode BVView_Default(BV bv,PetscViewer viewer)496: {
497: PetscErrorCode ierr;
498: PetscInt j;
499: Vec v;
500: PetscViewerFormat format;
501: PetscBool isascii,ismatlab=PETSC_FALSE;
502: const char *bvname,*name;
505: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
506: if (isascii) {
507: PetscViewerGetFormat(viewer,&format);
508: if (format == PETSC_VIEWER_ASCII_MATLAB) ismatlab = PETSC_TRUE;
509: }
510: if (ismatlab) {
511: PetscObjectGetName((PetscObject)bv,&bvname);
512: PetscViewerASCIIPrintf(viewer,"%s=[];\n",bvname);
513: }
514: for (j=bv->nc;j<bv->nc+bv->m;j++) {
515: BVGetColumn(bv,j,&v);
516: VecView(v,viewer);
517: if (ismatlab) {
518: PetscObjectGetName((PetscObject)v,&name);
519: PetscViewerASCIIPrintf(viewer,"%s=[%s,%s];clear %s\n",bvname,bvname,name,name);
520: }
521: BVRestoreColumn(bv,j,&v);
522: }
523: return(0);
524: }
528: /*@C
529: BVView - Prints the BV data structure.
531: Collective on BV533: Input Parameters:
534: + bv - the BV context
535: - viewer - optional visualization context
537: Note:
538: The available visualization contexts include
539: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
540: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
541: output where only the first processor opens
542: the file. All other processors send their
543: data to the first processor to print.
545: The user can open an alternative visualization contexts with
546: PetscViewerASCIIOpen() (output to a specified file).
548: Level: beginner
550: .seealso: PetscViewerASCIIOpen()
551: @*/
552: PetscErrorCode BVView(BV bv,PetscViewer viewer)553: {
554: PetscErrorCode ierr;
555: PetscBool isascii;
556: PetscViewerFormat format;
557: const char *orthname[2] = {"classical","modified"};
558: const char *refname[3] = {"if needed","never","always"};
559: const char *borthname[2] = {"Gram-Schmidt","Cholesky"};
563: if (!viewer) {
564: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)bv),&viewer);
565: }
568: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
569: if (isascii) {
570: PetscObjectPrintClassNamePrefixType((PetscObject)bv,viewer);
571: PetscViewerGetFormat(viewer,&format);
572: PetscViewerASCIIPushTab(viewer);
573: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
574: PetscViewerASCIIPrintf(viewer,"%D columns of global length %D\n",bv->m,bv->N);
575: if (bv->nc>0) {
576: PetscViewerASCIIPrintf(viewer,"number of constraints: %D\n",bv->nc);
577: }
578: PetscViewerASCIIPrintf(viewer,"vector orthogonalization method: %s Gram-Schmidt\n",orthname[bv->orthog_type]);
579: switch (bv->orthog_ref) {
580: case BV_ORTHOG_REFINE_IFNEEDED:
581: PetscViewerASCIIPrintf(viewer,"orthogonalization refinement: %s (eta: %g)\n",refname[bv->orthog_ref],(double)bv->orthog_eta);
582: break;
583: case BV_ORTHOG_REFINE_NEVER:
584: case BV_ORTHOG_REFINE_ALWAYS:
585: PetscViewerASCIIPrintf(viewer,"orthogonalization refinement: %s\n",refname[bv->orthog_ref]);
586: break;
587: }
588: PetscViewerASCIIPrintf(viewer,"block orthogonalization method: %s\n",borthname[bv->orthog_block]);
589: if (bv->matrix) {
590: if (bv->indef) {
591: PetscViewerASCIIPrintf(viewer,"indefinite inner product\n");
592: } else {
593: PetscViewerASCIIPrintf(viewer,"non-standard inner product\n");
594: }
595: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
596: MatView(bv->matrix,viewer);
597: PetscViewerPopFormat(viewer);
598: }
599: switch (bv->vmm) {
600: case BV_MATMULT_VECS:
601: PetscViewerASCIIPrintf(viewer,"doing matmult as matrix-vector products\n");
602: break;
603: case BV_MATMULT_MAT:
604: PetscViewerASCIIPrintf(viewer,"doing matmult as a single matrix-matrix product\n");
605: break;
606: case BV_MATMULT_MAT_SAVE:
607: PetscViewerASCIIPrintf(viewer,"doing matmult as a single matrix-matrix product, saving aux matrices\n");
608: break;
609: }
610: } else {
611: if (bv->ops->view) { (*bv->ops->view)(bv,viewer); }
612: else { BVView_Default(bv,viewer); }
613: }
614: PetscViewerASCIIPopTab(viewer);
615: } else {
616: (*bv->ops->view)(bv,viewer);
617: }
618: return(0);
619: }
623: /*@C
624: BVRegister - Adds a new storage format to de BV package.
626: Not collective
628: Input Parameters:
629: + name - name of a new user-defined BV630: - function - routine to create context
632: Notes:
633: BVRegister() may be called multiple times to add several user-defined
634: basis vectors.
636: Level: advanced
638: .seealso: BVRegisterAll()
639: @*/
640: PetscErrorCode BVRegister(const char *name,PetscErrorCode (*function)(BV))641: {
645: PetscFunctionListAdd(&BVList,name,function);
646: return(0);
647: }
651: PetscErrorCode BVAllocateWork_Private(BV bv,PetscInt s)652: {
656: if (s>bv->lwork) {
657: PetscFree(bv->work);
658: PetscMalloc1(s,&bv->work);
659: PetscLogObjectMemory((PetscObject)bv,(s-bv->lwork)*sizeof(PetscScalar));
660: bv->lwork = s;
661: }
662: return(0);
663: }