Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 : /*
11 : BV (basis vectors) interface routines, callable by users
12 : */
13 :
14 : #include <slepc/private/bvimpl.h> /*I "slepcbv.h" I*/
15 :
16 : PetscClassId BV_CLASSID = 0;
17 : 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_Normalize = 0,BV_SetRandom = 0,BV_MatMult = 0,BV_MatMultVec = 0,BV_MatProject = 0,BV_SVDAndRank = 0;
18 : static PetscBool BVPackageInitialized = PETSC_FALSE;
19 : MPI_Op MPIU_TSQR = 0,MPIU_LAPY2;
20 :
21 : const char *BVOrthogTypes[] = {"CGS","MGS","BVOrthogType","BV_ORTHOG_",NULL};
22 : const char *BVOrthogRefineTypes[] = {"IFNEEDED","NEVER","ALWAYS","BVOrthogRefineType","BV_ORTHOG_REFINE_",NULL};
23 : const char *BVOrthogBlockTypes[] = {"GS","CHOL","TSQR","TSQRCHOL","SVQB","BVOrthogBlockType","BV_ORTHOG_BLOCK_",NULL};
24 : const char *BVMatMultTypes[] = {"VECS","MAT","MAT_SAVE","BVMatMultType","BV_MATMULT_",NULL};
25 : const char *BVSVDMethods[] = {"REFINE","QR","QR_CAA","BVSVDMethod","BV_SVD_METHOD_",NULL};
26 :
27 : /*@C
28 : BVFinalizePackage - This function destroys everything in the Slepc interface
29 : to the BV package. It is called from SlepcFinalize().
30 :
31 : Level: developer
32 :
33 : .seealso: SlepcFinalize()
34 : @*/
35 1357 : PetscErrorCode BVFinalizePackage(void)
36 : {
37 1357 : PetscFunctionBegin;
38 1357 : PetscCall(PetscFunctionListDestroy(&BVList));
39 1357 : PetscCallMPI(MPI_Op_free(&MPIU_TSQR));
40 1357 : PetscCallMPI(MPI_Op_free(&MPIU_LAPY2));
41 1357 : BVPackageInitialized = PETSC_FALSE;
42 1357 : BVRegisterAllCalled = PETSC_FALSE;
43 1357 : PetscFunctionReturn(PETSC_SUCCESS);
44 : }
45 :
46 : /*@C
47 : BVInitializePackage - This function initializes everything in the BV package.
48 : It is called from PetscDLLibraryRegister() when using dynamic libraries, and
49 : on the first call to BVCreate() when using static libraries.
50 :
51 : Level: developer
52 :
53 : .seealso: SlepcInitialize()
54 : @*/
55 21162 : PetscErrorCode BVInitializePackage(void)
56 : {
57 21162 : char logList[256];
58 21162 : PetscBool opt,pkg;
59 21162 : PetscClassId classids[1];
60 :
61 21162 : PetscFunctionBegin;
62 21162 : if (BVPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
63 1357 : BVPackageInitialized = PETSC_TRUE;
64 : /* Register Classes */
65 1357 : PetscCall(PetscClassIdRegister("Basis Vectors",&BV_CLASSID));
66 : /* Register Constructors */
67 1357 : PetscCall(BVRegisterAll());
68 : /* Register Events */
69 1357 : PetscCall(PetscLogEventRegister("BVCreate",BV_CLASSID,&BV_Create));
70 1357 : PetscCall(PetscLogEventRegister("BVCopy",BV_CLASSID,&BV_Copy));
71 1357 : PetscCall(PetscLogEventRegister("BVMult",BV_CLASSID,&BV_Mult));
72 1357 : PetscCall(PetscLogEventRegister("BVMultVec",BV_CLASSID,&BV_MultVec));
73 1357 : PetscCall(PetscLogEventRegister("BVMultInPlace",BV_CLASSID,&BV_MultInPlace));
74 1357 : PetscCall(PetscLogEventRegister("BVDot",BV_CLASSID,&BV_Dot));
75 1357 : PetscCall(PetscLogEventRegister("BVDotVec",BV_CLASSID,&BV_DotVec));
76 1357 : PetscCall(PetscLogEventRegister("BVOrthogonalize",BV_CLASSID,&BV_Orthogonalize));
77 1357 : PetscCall(PetscLogEventRegister("BVOrthogonalizeV",BV_CLASSID,&BV_OrthogonalizeVec));
78 1357 : PetscCall(PetscLogEventRegister("BVScale",BV_CLASSID,&BV_Scale));
79 1357 : PetscCall(PetscLogEventRegister("BVNorm",BV_CLASSID,&BV_Norm));
80 1357 : PetscCall(PetscLogEventRegister("BVNormVec",BV_CLASSID,&BV_NormVec));
81 1357 : PetscCall(PetscLogEventRegister("BVNormalize",BV_CLASSID,&BV_Normalize));
82 1357 : PetscCall(PetscLogEventRegister("BVSetRandom",BV_CLASSID,&BV_SetRandom));
83 1357 : PetscCall(PetscLogEventRegister("BVMatMult",BV_CLASSID,&BV_MatMult));
84 1357 : PetscCall(PetscLogEventRegister("BVMatMultVec",BV_CLASSID,&BV_MatMultVec));
85 1357 : PetscCall(PetscLogEventRegister("BVMatProject",BV_CLASSID,&BV_MatProject));
86 1357 : PetscCall(PetscLogEventRegister("BVSVDAndRank",BV_CLASSID,&BV_SVDAndRank));
87 : /* MPI reduction operation used in BVOrthogonalize */
88 1357 : PetscCallMPI(MPI_Op_create(SlepcGivensPacked,PETSC_FALSE,&MPIU_TSQR));
89 1357 : PetscCallMPI(MPI_Op_create(SlepcPythag,PETSC_TRUE,&MPIU_LAPY2));
90 : /* Process Info */
91 1357 : classids[0] = BV_CLASSID;
92 1357 : PetscCall(PetscInfoProcessClass("bv",1,&classids[0]));
93 : /* Process summary exclusions */
94 1357 : PetscCall(PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt));
95 1357 : if (opt) {
96 10 : PetscCall(PetscStrInList("bv",logList,',',&pkg));
97 10 : if (pkg) PetscCall(PetscLogEventDeactivateClass(BV_CLASSID));
98 : }
99 : /* Register package finalizer */
100 1357 : PetscCall(PetscRegisterFinalize(BVFinalizePackage));
101 1357 : PetscFunctionReturn(PETSC_SUCCESS);
102 : }
103 :
104 : /*@
105 : BVDestroy - Destroys BV context that was created with BVCreate().
106 :
107 : Collective
108 :
109 : Input Parameter:
110 : . bv - the basis vectors context
111 :
112 : Level: beginner
113 :
114 : .seealso: BVCreate()
115 : @*/
116 59048 : PetscErrorCode BVDestroy(BV *bv)
117 : {
118 59048 : PetscFunctionBegin;
119 59048 : if (!*bv) PetscFunctionReturn(PETSC_SUCCESS);
120 14615 : PetscValidHeaderSpecific(*bv,BV_CLASSID,1);
121 14615 : PetscCheck(!(*bv)->lsplit,PetscObjectComm((PetscObject)*bv),PETSC_ERR_ARG_WRONGSTATE,"Must call BVRestoreSplit before destroying the BV");
122 14615 : if (--((PetscObject)*bv)->refct > 0) { *bv = NULL; PetscFunctionReturn(PETSC_SUCCESS); }
123 14376 : PetscTryTypeMethod(*bv,destroy);
124 14376 : PetscCall(PetscLayoutDestroy(&(*bv)->map));
125 14376 : PetscCall(PetscFree((*bv)->vtype));
126 14376 : PetscCall(MatDestroy(&(*bv)->matrix));
127 14376 : PetscCall(VecDestroy(&(*bv)->Bx));
128 14376 : PetscCall(VecDestroy(&(*bv)->buffer));
129 14376 : PetscCall(BVDestroy(&(*bv)->cached));
130 14376 : PetscCall(BVDestroy(&(*bv)->L));
131 14376 : PetscCall(BVDestroy(&(*bv)->R));
132 14376 : PetscCall(PetscFree((*bv)->work));
133 14376 : PetscCall(PetscFree2((*bv)->h,(*bv)->c));
134 14376 : PetscCall(VecDestroy(&(*bv)->omega));
135 14376 : PetscCall(MatDestroy(&(*bv)->Acreate));
136 14376 : PetscCall(MatDestroy(&(*bv)->Aget));
137 14376 : PetscCall(MatDestroy(&(*bv)->Abuffer));
138 14376 : PetscCall(PetscRandomDestroy(&(*bv)->rand));
139 14376 : PetscCall(PetscHeaderDestroy(bv));
140 14376 : PetscFunctionReturn(PETSC_SUCCESS);
141 : }
142 :
143 : /*@
144 : BVCreate - Creates a basis vectors context.
145 :
146 : Collective
147 :
148 : Input Parameter:
149 : . comm - MPI communicator
150 :
151 : Output Parameter:
152 : . newbv - location to put the basis vectors context
153 :
154 : Level: beginner
155 :
156 : .seealso: BVSetUp(), BVDestroy(), BV
157 : @*/
158 14376 : PetscErrorCode BVCreate(MPI_Comm comm,BV *newbv)
159 : {
160 14376 : BV bv;
161 :
162 14376 : PetscFunctionBegin;
163 14376 : PetscAssertPointer(newbv,2);
164 14376 : PetscCall(BVInitializePackage());
165 14376 : PetscCall(SlepcHeaderCreate(bv,BV_CLASSID,"BV","Basis Vectors","BV",comm,BVDestroy,BVView));
166 :
167 14376 : bv->map = NULL;
168 14376 : bv->vtype = NULL;
169 14376 : bv->n = -1;
170 14376 : bv->N = -1;
171 14376 : bv->m = 0;
172 14376 : bv->l = 0;
173 14376 : bv->k = 0;
174 14376 : bv->nc = 0;
175 14376 : bv->ld = 0;
176 14376 : bv->orthog_type = BV_ORTHOG_CGS;
177 14376 : bv->orthog_ref = BV_ORTHOG_REFINE_IFNEEDED;
178 14376 : bv->orthog_eta = 0.7071;
179 14376 : bv->orthog_block = BV_ORTHOG_BLOCK_GS;
180 14376 : bv->matrix = NULL;
181 14376 : bv->indef = PETSC_FALSE;
182 14376 : bv->vmm = BV_MATMULT_MAT;
183 14376 : bv->rrandom = PETSC_FALSE;
184 14376 : bv->deftol = 10*PETSC_MACHINE_EPSILON;
185 :
186 14376 : bv->buffer = NULL;
187 14376 : bv->Abuffer = NULL;
188 14376 : bv->Bx = NULL;
189 14376 : bv->xid = 0;
190 14376 : bv->xstate = 0;
191 14376 : bv->cv[0] = NULL;
192 14376 : bv->cv[1] = NULL;
193 14376 : bv->ci[0] = -1;
194 14376 : bv->ci[1] = -1;
195 14376 : bv->st[0] = -1;
196 14376 : bv->st[1] = -1;
197 14376 : bv->id[0] = 0;
198 14376 : bv->id[1] = 0;
199 14376 : bv->h = NULL;
200 14376 : bv->c = NULL;
201 14376 : bv->omega = NULL;
202 14376 : bv->defersfo = PETSC_FALSE;
203 14376 : bv->cached = NULL;
204 14376 : bv->bvstate = 0;
205 14376 : bv->L = NULL;
206 14376 : bv->R = NULL;
207 14376 : bv->lstate = 0;
208 14376 : bv->rstate = 0;
209 14376 : bv->lsplit = 0;
210 14376 : bv->issplit = 0;
211 14376 : bv->splitparent = NULL;
212 14376 : bv->rand = NULL;
213 14376 : bv->Acreate = NULL;
214 14376 : bv->Aget = NULL;
215 14376 : bv->cuda = PETSC_FALSE;
216 14376 : bv->hip = PETSC_FALSE;
217 14376 : bv->sfocalled = PETSC_FALSE;
218 14376 : bv->work = NULL;
219 14376 : bv->lwork = 0;
220 14376 : bv->data = NULL;
221 :
222 14376 : *newbv = bv;
223 14376 : PetscFunctionReturn(PETSC_SUCCESS);
224 : }
225 :
226 : /*@
227 : BVCreateFromMat - Creates a basis vectors object from a dense Mat object.
228 :
229 : Collective
230 :
231 : Input Parameter:
232 : . A - a dense tall-skinny matrix
233 :
234 : Output Parameter:
235 : . bv - the new basis vectors context
236 :
237 : Notes:
238 : The matrix values are copied to the BV data storage, memory is not shared.
239 :
240 : The communicator of the BV object will be the same as A, and so will be
241 : the dimensions.
242 :
243 : Level: intermediate
244 :
245 : .seealso: BVCreate(), BVDestroy(), BVCreateMat()
246 : @*/
247 74 : PetscErrorCode BVCreateFromMat(Mat A,BV *bv)
248 : {
249 74 : PetscInt n,N,k;
250 74 : VecType vtype;
251 :
252 74 : PetscFunctionBegin;
253 74 : PetscValidHeaderSpecific(A,MAT_CLASSID,1);
254 :
255 74 : PetscCall(MatGetSize(A,&N,&k));
256 74 : PetscCall(MatGetLocalSize(A,&n,NULL));
257 74 : PetscCall(MatGetVecType(A,&vtype));
258 74 : PetscCall(BVCreate(PetscObjectComm((PetscObject)A),bv));
259 74 : PetscCall(BVSetSizes(*bv,n,N,k));
260 74 : PetscCall(BVSetVecType(*bv,vtype));
261 :
262 74 : (*bv)->Acreate = A;
263 74 : PetscCall(PetscObjectReference((PetscObject)A));
264 74 : PetscFunctionReturn(PETSC_SUCCESS);
265 : }
266 :
267 : /*@
268 : BVInsertVec - Insert a vector into the specified column.
269 :
270 : Logically Collective
271 :
272 : Input Parameters:
273 : + V - basis vectors
274 : . j - the column of V to be overwritten
275 : - w - the vector to be copied
276 :
277 : Level: intermediate
278 :
279 : .seealso: BVInsertVecs()
280 : @*/
281 54739 : PetscErrorCode BVInsertVec(BV V,PetscInt j,Vec w)
282 : {
283 54739 : PetscInt n,N;
284 54739 : Vec v;
285 :
286 54739 : PetscFunctionBegin;
287 54739 : PetscValidHeaderSpecific(V,BV_CLASSID,1);
288 164217 : PetscValidLogicalCollectiveInt(V,j,2);
289 54739 : PetscValidHeaderSpecific(w,VEC_CLASSID,3);
290 54739 : PetscValidType(V,1);
291 54739 : BVCheckSizes(V,1);
292 54739 : PetscCheckSameComm(V,1,w,3);
293 :
294 54739 : PetscCall(VecGetSize(w,&N));
295 54739 : PetscCall(VecGetLocalSize(w,&n));
296 54739 : PetscCheck(N==V->N && n==V->n,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Vec sizes (global %" PetscInt_FMT ", local %" PetscInt_FMT ") do not match BV sizes (global %" PetscInt_FMT ", local %" PetscInt_FMT ")",N,n,V->N,V->n);
297 54739 : PetscCheck(j>=-V->nc && j<V->m,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument j has wrong value %" PetscInt_FMT ", should be between %" PetscInt_FMT " and %" PetscInt_FMT,j,-V->nc,V->m-1);
298 :
299 54739 : PetscCall(BVGetColumn(V,j,&v));
300 54739 : PetscCall(VecCopy(w,v));
301 54739 : PetscCall(BVRestoreColumn(V,j,&v));
302 54739 : PetscCall(PetscObjectStateIncrease((PetscObject)V));
303 54739 : PetscFunctionReturn(PETSC_SUCCESS);
304 : }
305 :
306 : /*@
307 : BVInsertVecs - Insert a set of vectors into the specified columns.
308 :
309 : Collective
310 :
311 : Input Parameters:
312 : + V - basis vectors
313 : . s - first column of V to be overwritten
314 : . W - set of vectors to be copied
315 : - orth - flag indicating if the vectors must be orthogonalized
316 :
317 : Input/Output Parameter:
318 : . m - number of input vectors, on output the number of linearly independent
319 : vectors
320 :
321 : Notes:
322 : Copies the contents of vectors W to V(:,s:s+n). If the orthogonalization
323 : flag is set, then the vectors are copied one by one and then orthogonalized
324 : against the previous ones. If any of them is linearly dependent then it
325 : is discarded and the value of m is decreased.
326 :
327 : Level: intermediate
328 :
329 : .seealso: BVInsertVec(), BVOrthogonalizeColumn()
330 : @*/
331 328 : PetscErrorCode BVInsertVecs(BV V,PetscInt s,PetscInt *m,Vec *W,PetscBool orth)
332 : {
333 328 : PetscInt n,N,i,ndep;
334 328 : PetscBool lindep;
335 328 : PetscReal norm;
336 328 : Vec v;
337 :
338 328 : PetscFunctionBegin;
339 328 : PetscValidHeaderSpecific(V,BV_CLASSID,1);
340 984 : PetscValidLogicalCollectiveInt(V,s,2);
341 328 : PetscAssertPointer(m,3);
342 984 : PetscValidLogicalCollectiveInt(V,*m,3);
343 328 : if (!*m) PetscFunctionReturn(PETSC_SUCCESS);
344 328 : PetscCheck(*m>0,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %" PetscInt_FMT ") cannot be negative",*m);
345 328 : PetscAssertPointer(W,4);
346 328 : PetscValidHeaderSpecific(*W,VEC_CLASSID,4);
347 984 : PetscValidLogicalCollectiveBool(V,orth,5);
348 328 : PetscValidType(V,1);
349 328 : BVCheckSizes(V,1);
350 328 : PetscCheckSameComm(V,1,*W,4);
351 :
352 328 : PetscCall(VecGetSize(*W,&N));
353 328 : PetscCall(VecGetLocalSize(*W,&n));
354 328 : PetscCheck(N==V->N && n==V->n,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Vec sizes (global %" PetscInt_FMT ", local %" PetscInt_FMT ") do not match BV sizes (global %" PetscInt_FMT ", local %" PetscInt_FMT ")",N,n,V->N,V->n);
355 328 : PetscCheck(s>=0 && s<V->m,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument s has wrong value %" PetscInt_FMT ", should be between 0 and %" PetscInt_FMT,s,V->m-1);
356 328 : PetscCheck(s+(*m)<=V->m,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Too many vectors provided, there is only room for %" PetscInt_FMT,V->m);
357 :
358 : ndep = 0;
359 722 : for (i=0;i<*m;i++) {
360 394 : PetscCall(BVGetColumn(V,s+i-ndep,&v));
361 394 : PetscCall(VecCopy(W[i],v));
362 394 : PetscCall(BVRestoreColumn(V,s+i-ndep,&v));
363 394 : if (orth) {
364 394 : PetscCall(BVOrthogonalizeColumn(V,s+i-ndep,NULL,&norm,&lindep));
365 394 : if (norm==0.0 || lindep) {
366 0 : PetscCall(PetscInfo(V,"Removing linearly dependent vector %" PetscInt_FMT "\n",i));
367 0 : ndep++;
368 394 : } else PetscCall(BVScaleColumn(V,s+i-ndep,1.0/norm));
369 : }
370 : }
371 328 : *m -= ndep;
372 328 : PetscCall(PetscObjectStateIncrease((PetscObject)V));
373 328 : PetscFunctionReturn(PETSC_SUCCESS);
374 : }
375 :
376 : /*@
377 : BVInsertConstraints - Insert a set of vectors as constraints.
378 :
379 : Collective
380 :
381 : Input Parameters:
382 : + V - basis vectors
383 : - C - set of vectors to be inserted as constraints
384 :
385 : Input/Output Parameter:
386 : . nc - number of input vectors, on output the number of linearly independent
387 : vectors
388 :
389 : Notes:
390 : The constraints are relevant only during orthogonalization. Constraint
391 : vectors span a subspace that is deflated in every orthogonalization
392 : operation, so they are intended for removing those directions from the
393 : orthogonal basis computed in regular BV columns.
394 :
395 : Constraints are not stored in regular BV columns, but in a special part of
396 : the storage. They can be accessed with negative indices in BVGetColumn().
397 :
398 : This operation is DESTRUCTIVE, meaning that all data contained in the
399 : columns of V is lost. This is typically invoked just after creating the BV.
400 : Once a set of constraints has been set, it is not allowed to call this
401 : function again.
402 :
403 : The vectors are copied one by one and then orthogonalized against the
404 : previous ones. If any of them is linearly dependent then it is discarded
405 : and the value of nc is decreased. The behaviour is similar to BVInsertVecs().
406 :
407 : Level: advanced
408 :
409 : .seealso: BVInsertVecs(), BVOrthogonalizeColumn(), BVGetColumn(), BVGetNumConstraints()
410 : @*/
411 44 : PetscErrorCode BVInsertConstraints(BV V,PetscInt *nc,Vec *C)
412 : {
413 44 : PetscInt msave;
414 :
415 44 : PetscFunctionBegin;
416 44 : PetscValidHeaderSpecific(V,BV_CLASSID,1);
417 44 : PetscAssertPointer(nc,2);
418 132 : PetscValidLogicalCollectiveInt(V,*nc,2);
419 44 : if (!*nc) PetscFunctionReturn(PETSC_SUCCESS);
420 44 : PetscCheck(*nc>0,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Number of constraints (given %" PetscInt_FMT ") cannot be negative",*nc);
421 44 : PetscAssertPointer(C,3);
422 44 : PetscValidHeaderSpecific(*C,VEC_CLASSID,3);
423 44 : PetscValidType(V,1);
424 44 : BVCheckSizes(V,1);
425 44 : PetscCheckSameComm(V,1,*C,3);
426 44 : PetscCheck(!V->issplit,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_WRONGSTATE,"Operation not permitted for a BV obtained from BVGetSplit");
427 44 : PetscCheck(!V->nc,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_WRONGSTATE,"Constraints already present in this BV object");
428 44 : PetscCheck(V->ci[0]==-1 && V->ci[1]==-1,PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Cannot call BVInsertConstraints after BVGetColumn");
429 :
430 44 : msave = V->m;
431 44 : PetscCall(BVResize(V,*nc+V->m,PETSC_FALSE));
432 44 : PetscCall(BVInsertVecs(V,0,nc,C,PETSC_TRUE));
433 44 : V->nc = *nc;
434 44 : V->m = msave;
435 44 : V->ci[0] = -V->nc-1;
436 44 : V->ci[1] = -V->nc-1;
437 44 : PetscCall(PetscObjectStateIncrease((PetscObject)V));
438 44 : PetscFunctionReturn(PETSC_SUCCESS);
439 : }
440 :
441 : /*@
442 : BVSetOptionsPrefix - Sets the prefix used for searching for all
443 : BV options in the database.
444 :
445 : Logically Collective
446 :
447 : Input Parameters:
448 : + bv - the basis vectors context
449 : - prefix - the prefix string to prepend to all BV option requests
450 :
451 : Notes:
452 : A hyphen (-) must NOT be given at the beginning of the prefix name.
453 : The first character of all runtime options is AUTOMATICALLY the
454 : hyphen.
455 :
456 : Level: advanced
457 :
458 : .seealso: BVAppendOptionsPrefix(), BVGetOptionsPrefix()
459 : @*/
460 241 : PetscErrorCode BVSetOptionsPrefix(BV bv,const char *prefix)
461 : {
462 241 : PetscFunctionBegin;
463 241 : PetscValidHeaderSpecific(bv,BV_CLASSID,1);
464 241 : PetscCall(PetscObjectSetOptionsPrefix((PetscObject)bv,prefix));
465 241 : PetscFunctionReturn(PETSC_SUCCESS);
466 : }
467 :
468 : /*@
469 : BVAppendOptionsPrefix - Appends to the prefix used for searching for all
470 : BV options in the database.
471 :
472 : Logically Collective
473 :
474 : Input Parameters:
475 : + bv - the basis vectors context
476 : - prefix - the prefix string to prepend to all BV option requests
477 :
478 : Notes:
479 : A hyphen (-) must NOT be given at the beginning of the prefix name.
480 : The first character of all runtime options is AUTOMATICALLY the
481 : hyphen.
482 :
483 : Level: advanced
484 :
485 : .seealso: BVSetOptionsPrefix(), BVGetOptionsPrefix()
486 : @*/
487 201 : PetscErrorCode BVAppendOptionsPrefix(BV bv,const char *prefix)
488 : {
489 201 : PetscFunctionBegin;
490 201 : PetscValidHeaderSpecific(bv,BV_CLASSID,1);
491 201 : PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)bv,prefix));
492 201 : PetscFunctionReturn(PETSC_SUCCESS);
493 : }
494 :
495 : /*@
496 : BVGetOptionsPrefix - Gets the prefix used for searching for all
497 : BV options in the database.
498 :
499 : Not Collective
500 :
501 : Input Parameters:
502 : . bv - the basis vectors context
503 :
504 : Output Parameters:
505 : . prefix - pointer to the prefix string used, is returned
506 :
507 : Note:
508 : On the Fortran side, the user should pass in a string 'prefix' of
509 : sufficient length to hold the prefix.
510 :
511 : Level: advanced
512 :
513 : .seealso: BVSetOptionsPrefix(), BVAppendOptionsPrefix()
514 : @*/
515 0 : PetscErrorCode BVGetOptionsPrefix(BV bv,const char *prefix[])
516 : {
517 0 : PetscFunctionBegin;
518 0 : PetscValidHeaderSpecific(bv,BV_CLASSID,1);
519 0 : PetscAssertPointer(prefix,2);
520 0 : PetscCall(PetscObjectGetOptionsPrefix((PetscObject)bv,prefix));
521 0 : PetscFunctionReturn(PETSC_SUCCESS);
522 : }
523 :
524 : /*@
525 : BVView - Prints the BV data structure.
526 :
527 : Collective
528 :
529 : Input Parameters:
530 : + bv - the BV context
531 : - viewer - optional visualization context
532 :
533 : Note:
534 : The available visualization contexts include
535 : + PETSC_VIEWER_STDOUT_SELF - standard output (default)
536 : - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
537 : output where only the first processor opens
538 : the file. All other processors send their
539 : data to the first processor to print.
540 :
541 : The user can open an alternative visualization contexts with
542 : PetscViewerASCIIOpen() (output to a specified file).
543 :
544 : Level: beginner
545 :
546 : .seealso: BVCreate()
547 : @*/
548 74 : PetscErrorCode BVView(BV bv,PetscViewer viewer)
549 : {
550 74 : PetscBool isascii;
551 74 : PetscViewerFormat format;
552 74 : const char *orthname[2] = {"classical","modified"};
553 74 : const char *refname[3] = {"if needed","never","always"};
554 :
555 74 : PetscFunctionBegin;
556 74 : PetscValidHeaderSpecific(bv,BV_CLASSID,1);
557 74 : if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)bv),&viewer));
558 74 : PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
559 :
560 74 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
561 74 : if (isascii) {
562 74 : PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)bv,viewer));
563 74 : PetscCall(PetscViewerGetFormat(viewer,&format));
564 74 : if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
565 68 : PetscCall(PetscViewerASCIIPrintf(viewer," %" PetscInt_FMT " columns of global length %" PetscInt_FMT "%s\n",bv->m,bv->N,bv->cuda?" (CUDA)":bv->hip?" (HIP)":""));
566 34 : if (bv->nc>0) PetscCall(PetscViewerASCIIPrintf(viewer," number of constraints: %" PetscInt_FMT "\n",bv->nc));
567 34 : PetscCall(PetscViewerASCIIPrintf(viewer," vector orthogonalization method: %s Gram-Schmidt\n",orthname[bv->orthog_type]));
568 34 : switch (bv->orthog_ref) {
569 34 : case BV_ORTHOG_REFINE_IFNEEDED:
570 34 : PetscCall(PetscViewerASCIIPrintf(viewer," orthogonalization refinement: %s (eta: %g)\n",refname[bv->orthog_ref],(double)bv->orthog_eta));
571 : break;
572 0 : case BV_ORTHOG_REFINE_NEVER:
573 : case BV_ORTHOG_REFINE_ALWAYS:
574 0 : PetscCall(PetscViewerASCIIPrintf(viewer," orthogonalization refinement: %s\n",refname[bv->orthog_ref]));
575 : break;
576 : }
577 34 : PetscCall(PetscViewerASCIIPrintf(viewer," block orthogonalization method: %s\n",BVOrthogBlockTypes[bv->orthog_block]));
578 34 : if (bv->matrix) {
579 0 : if (bv->indef) PetscCall(PetscViewerASCIIPrintf(viewer," indefinite inner product\n"));
580 0 : else PetscCall(PetscViewerASCIIPrintf(viewer," non-standard inner product\n"));
581 0 : PetscCall(PetscViewerASCIIPrintf(viewer," tolerance for definite inner product: %g\n",(double)bv->deftol));
582 0 : PetscCall(PetscViewerASCIIPrintf(viewer," inner product matrix:\n"));
583 0 : PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO));
584 0 : PetscCall(PetscViewerASCIIPushTab(viewer));
585 0 : PetscCall(MatView(bv->matrix,viewer));
586 0 : PetscCall(PetscViewerASCIIPopTab(viewer));
587 0 : PetscCall(PetscViewerPopFormat(viewer));
588 : }
589 34 : switch (bv->vmm) {
590 4 : case BV_MATMULT_VECS:
591 4 : PetscCall(PetscViewerASCIIPrintf(viewer," doing matmult as matrix-vector products\n"));
592 : break;
593 30 : case BV_MATMULT_MAT:
594 30 : PetscCall(PetscViewerASCIIPrintf(viewer," doing matmult as a single matrix-matrix product\n"));
595 : break;
596 0 : case BV_MATMULT_MAT_SAVE:
597 0 : PetscCall(PetscViewerASCIIPrintf(viewer," mat_save is deprecated, use mat\n"));
598 : break;
599 : }
600 34 : if (bv->rrandom) PetscCall(PetscViewerASCIIPrintf(viewer," generating random vectors independent of the number of processes\n"));
601 : }
602 : }
603 74 : PetscTryTypeMethod(bv,view,viewer);
604 74 : PetscFunctionReturn(PETSC_SUCCESS);
605 : }
606 :
607 : /*@
608 : BVViewFromOptions - View from options
609 :
610 : Collective
611 :
612 : Input Parameters:
613 : + bv - the basis vectors context
614 : . obj - optional object
615 : - name - command line option
616 :
617 : Level: intermediate
618 :
619 : .seealso: BVView(), BVCreate()
620 : @*/
621 0 : PetscErrorCode BVViewFromOptions(BV bv,PetscObject obj,const char name[])
622 : {
623 0 : PetscFunctionBegin;
624 0 : PetscValidHeaderSpecific(bv,BV_CLASSID,1);
625 0 : PetscCall(PetscObjectViewFromOptions((PetscObject)bv,obj,name));
626 0 : PetscFunctionReturn(PETSC_SUCCESS);
627 : }
628 :
629 : /*@C
630 : BVRegister - Adds a new storage format to the BV package.
631 :
632 : Not Collective
633 :
634 : Input Parameters:
635 : + name - name of a new user-defined BV
636 : - function - routine to create context
637 :
638 : Notes:
639 : BVRegister() may be called multiple times to add several user-defined
640 : basis vectors.
641 :
642 : Level: advanced
643 :
644 : .seealso: BVRegisterAll()
645 : @*/
646 6785 : PetscErrorCode BVRegister(const char *name,PetscErrorCode (*function)(BV))
647 : {
648 6785 : PetscFunctionBegin;
649 6785 : PetscCall(BVInitializePackage());
650 6785 : PetscCall(PetscFunctionListAdd(&BVList,name,function));
651 6785 : PetscFunctionReturn(PETSC_SUCCESS);
652 : }
653 :
654 61155 : PetscErrorCode BVAllocateWork_Private(BV bv,PetscInt s)
655 : {
656 61155 : PetscFunctionBegin;
657 61155 : if (s>bv->lwork) {
658 15611 : PetscCall(PetscFree(bv->work));
659 15611 : PetscCall(PetscMalloc1(s,&bv->work));
660 15611 : bv->lwork = s;
661 : }
662 61155 : PetscFunctionReturn(PETSC_SUCCESS);
663 : }
664 :
665 : #if defined(PETSC_USE_DEBUG) && !defined(PETSC_CLANG_STATIC_ANALYZER)
666 : /*
667 : SlepcDebugBVView - partially view a BV object, to be used from within a debugger.
668 :
669 : ini, end: columns to be viewed
670 : s: name of Matlab variable
671 : filename: optionally write output to a file
672 : */
673 0 : PETSC_UNUSED PetscErrorCode SlepcDebugBVView(BV bv,PetscInt ini,PetscInt end,const char *s,const char *filename)
674 : {
675 0 : PetscInt N,m;
676 0 : PetscScalar *array;
677 :
678 0 : PetscFunctionBegin;
679 0 : PetscCall(BVGetArray(bv,&array));
680 0 : PetscCall(BVGetSizes(bv,NULL,&N,&m));
681 0 : PetscCall(SlepcDebugViewMatrix(N,end-ini+1,array+ini*N,NULL,bv->ld,s,filename));
682 0 : PetscCall(BVRestoreArray(bv,&array));
683 0 : PetscFunctionReturn(PETSC_SUCCESS);
684 : }
685 : #endif
|