1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: #include <slepc/private/vecimplslepc.h> 13: /*@
14: VecNormalizeComplex - Normalizes a possibly complex vector by the 2-norm.
16: Collective on xr
18: Input parameters:
19: + xr - the real part of the vector (overwritten on output)
20: . xi - the imaginary part of the vector (not referenced if iscomplex is false)
21: - iscomplex - a flag indicating if the vector is complex
23: Output parameter:
24: . norm - the vector norm before normalization (can be set to NULL)
26: Level: developer
27: @*/
28: PetscErrorCode VecNormalizeComplex(Vec xr,Vec xi,PetscBool iscomplex,PetscReal *norm) 29: {
31: #if !defined(PETSC_USE_COMPLEX)
32: PetscReal normr,normi,alpha;
33: #endif
37: #if !defined(PETSC_USE_COMPLEX)
38: if (iscomplex) {
40: VecNormBegin(xr,NORM_2,&normr);
41: VecNormBegin(xi,NORM_2,&normi);
42: VecNormEnd(xr,NORM_2,&normr);
43: VecNormEnd(xi,NORM_2,&normi);
44: alpha = SlepcAbsEigenvalue(normr,normi);
45: if (norm) *norm = alpha;
46: alpha = 1.0 / alpha;
47: VecScale(xr,alpha);
48: VecScale(xi,alpha);
49: } else
50: #endif
51: {
52: VecNormalize(xr,norm);
53: }
54: return(0);
55: }
57: static PetscErrorCode VecCheckOrthogonality_Private(Vec V[],PetscInt nv,Vec W[],PetscInt nw,Mat B,PetscViewer viewer,PetscReal *lev,PetscBool norm) 58: {
60: PetscInt i,j;
61: PetscScalar *vals;
62: PetscBool isascii;
63: Vec w;
66: if (!lev) {
67: if (!viewer) {
68: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)*V),&viewer);
69: }
72: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
73: if (!isascii) return(0);
74: }
76: PetscMalloc1(nv,&vals);
77: if (B) {
78: VecDuplicate(V[0],&w);
79: }
80: if (lev) *lev = 0.0;
81: for (i=0;i<nw;i++) {
82: if (B) {
83: if (W) {
84: MatMultTranspose(B,W[i],w);
85: } else {
86: MatMultTranspose(B,V[i],w);
87: }
88: } else {
89: if (W) w = W[i];
90: else w = V[i];
91: }
92: VecMDot(w,nv,V,vals);
93: for (j=0;j<nv;j++) {
94: if (lev) {
95: if (i!=j) *lev = PetscMax(*lev,PetscAbsScalar(vals[j]));
96: else if (norm) *lev = PetscMax(*lev,PetscAbsScalar(vals[j]-PetscRealConstant(1.0)));
97: } else {
98: #if !defined(PETSC_USE_COMPLEX)
99: PetscViewerASCIIPrintf(viewer," %12g ",(double)vals[j]);
100: #else
101: PetscViewerASCIIPrintf(viewer," %12g%+12gi ",(double)PetscRealPart(vals[j]),(double)PetscImaginaryPart(vals[j]));
102: #endif
103: }
104: }
105: if (!lev) { PetscViewerASCIIPrintf(viewer,"\n"); }
106: }
107: PetscFree(vals);
108: if (B) {
109: VecDestroy(&w);
110: }
111: return(0);
112: }
114: /*@
115: VecCheckOrthogonality - Checks (or prints) the level of (bi-)orthogonality
116: of a set of vectors.
118: Collective on V
120: Input parameters:
121: + V - a set of vectors
122: . nv - number of V vectors
123: . W - an alternative set of vectors (optional)
124: . nw - number of W vectors
125: . B - matrix defining the inner product (optional)
126: - viewer - optional visualization context
128: Output parameter:
129: . lev - level of orthogonality (optional)
131: Notes:
132: This function computes W'*V and prints the result. It is intended to check
133: the level of bi-orthogonality of the vectors in the two sets. If W is equal
134: to NULL then V is used, thus checking the orthogonality of the V vectors.
136: If matrix B is provided then the check uses the B-inner product, W'*B*V.
138: If lev is not NULL, it will contain the maximum entry of matrix
139: W'*V - I (in absolute value) omitting the diagonal. Otherwise, the matrix W'*V
140: is printed.
142: Level: developer
144: .seealso: VecCheckOrthonormality()
145: @*/
146: PetscErrorCode VecCheckOrthogonality(Vec V[],PetscInt nv,Vec W[],PetscInt nw,Mat B,PetscViewer viewer,PetscReal *lev)147: {
155: if (nv<=0 || nw<=0) return(0);
156: if (W) {
160: }
161: VecCheckOrthogonality_Private(V,nv,W,nw,B,viewer,lev,PETSC_FALSE);
162: return(0);
163: }
165: /*@
166: VecCheckOrthonormality - Checks (or prints) the level of (bi-)orthonormality
167: of a set of vectors.
169: Collective on V
171: Input parameters:
172: + V - a set of vectors
173: . nv - number of V vectors
174: . W - an alternative set of vectors (optional)
175: . nw - number of W vectors
176: . B - matrix defining the inner product (optional)
177: - viewer - optional visualization context
179: Output parameter:
180: . lev - level of orthogonality (optional)
182: Notes:
183: This function is equivalent to VecCheckOrthonormality(), but in addition it checks
184: that the diagonal of W'*V (or W'*B*V) is equal to all ones.
186: Level: developer
188: .seealso: VecCheckOrthogonality()
189: @*/
190: PetscErrorCode VecCheckOrthonormality(Vec V[],PetscInt nv,Vec W[],PetscInt nw,Mat B,PetscViewer viewer,PetscReal *lev)191: {
199: if (nv<=0 || nw<=0) return(0);
200: if (W) {
204: }
205: VecCheckOrthogonality_Private(V,nv,W,nw,B,viewer,lev,PETSC_TRUE);
206: return(0);
207: }
209: /*@
210: VecDuplicateEmpty - Creates a new vector of the same type as an existing vector,
211: but without internal array.
213: Collective on v
215: Input parameters:
216: . v - a vector to mimic
218: Output parameter:
219: . newv - location to put new vector
221: Note:
222: This is similar to VecDuplicate(), but the new vector does not have an internal
223: array, so the intended usage is with VecPlaceArray().
225: Level: developer
226: @*/
227: PetscErrorCode VecDuplicateEmpty(Vec v,Vec *newv)228: {
230: PetscBool standard,cuda,mpi;
231: PetscInt N,nloc,bs;
238: PetscObjectTypeCompareAny((PetscObject)v,&standard,VECSEQ,VECMPI,"");
239: PetscObjectTypeCompareAny((PetscObject)v,&cuda,VECSEQCUDA,VECMPICUDA,"");
240: if (standard || cuda) {
241: PetscObjectTypeCompareAny((PetscObject)v,&mpi,VECMPI,VECMPICUDA,"");
242: VecGetLocalSize(v,&nloc);
243: VecGetSize(v,&N);
244: VecGetBlockSize(v,&bs);
245: if (cuda) {
246: #if defined(PETSC_HAVE_CUDA)
247: if (mpi) {
248: VecCreateMPICUDAWithArray(PetscObjectComm((PetscObject)v),bs,nloc,N,NULL,newv);
249: } else {
250: VecCreateSeqCUDAWithArray(PetscObjectComm((PetscObject)v),bs,N,NULL,newv);
251: }
252: #endif
253: } else {
254: if (mpi) {
255: VecCreateMPIWithArray(PetscObjectComm((PetscObject)v),bs,nloc,N,NULL,newv);
256: } else {
257: VecCreateSeqWithArray(PetscObjectComm((PetscObject)v),bs,N,NULL,newv);
258: }
259: }
260: } else { /* standard duplicate, with internal array */
261: VecDuplicate(v,newv);
262: }
263: return(0);
264: }
266: /*@
267: VecSetRandomNormal - Sets all components of a vector to normally distributed random values.
269: Logically Collective on v
271: Input parameters:
272: + v - the vector to be filled with random values
273: . rctx - the random number context (can be NULL)
274: - w1,w2 - two work vectors (can be NULL)
276: Output parameter:
277: . v - the vector
279: Notes:
280: Fills the two work vectors with uniformly distributed random values (VecSetRandom)
281: and then applies the Box-Muller transform to get normally distributed values on v.
283: Level: developer
284: @*/
285: PetscErrorCode VecSetRandomNormal(Vec v,PetscRandom rctx,Vec w1,Vec w2)286: {
287: PetscErrorCode ierr;
288: const PetscScalar *x,*y;
289: PetscScalar *z;
290: PetscInt n,i;
291: PetscRandom rand=NULL;
292: Vec v1=NULL,v2=NULL;
301: if (!rctx) {
302: PetscRandomCreate(PetscObjectComm((PetscObject)v),&rand);
303: PetscRandomSetFromOptions(rand);
304: rctx = rand;
305: }
306: if (!w1) {
307: VecDuplicate(v,&v1);
308: w1 = v1;
309: }
310: if (!w2) {
311: VecDuplicate(v,&v2);
312: w2 = v2;
313: }
317: VecSetRandom(w1,rctx);
318: VecSetRandom(w2,rctx);
319: VecGetLocalSize(v,&n);
320: VecGetArrayWrite(v,&z);
321: VecGetArrayRead(w1,&x);
322: VecGetArrayRead(w2,&y);
323: for (i=0;i<n;i++) {
324: #if defined(PETSC_USE_COMPLEX)
325: z[i] = PetscCMPLX(PetscSqrtReal(-2.0*PetscLogReal(PetscRealPart(x[i])))*PetscCosReal(2.0*PETSC_PI*PetscRealPart(y[i])),PetscSqrtReal(-2.0*PetscLogReal(PetscImaginaryPart(x[i])))*PetscCosReal(2.0*PETSC_PI*PetscImaginaryPart(y[i])));
326: #else
327: z[i] = PetscSqrtReal(-2.0*PetscLogReal(x[i]))*PetscCosReal(2.0*PETSC_PI*y[i]);
328: #endif
329: }
330: VecRestoreArrayWrite(v,&z);
331: VecRestoreArrayRead(w1,&x);
332: VecRestoreArrayRead(w2,&y);
334: VecDestroy(&v1);
335: VecDestroy(&v2);
336: if (!rctx) {
337: PetscRandomDestroy(&rand);
338: }
339: return(0);
340: }