Actual source code: vecutil.c

slepc-3.8.2 2017-12-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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>       /*I "slepcvec.h" I*/

 13: /*@
 14:    VecNormalizeComplex - Normalizes a possibly complex vector by the 2-norm.

 16:    Collective on Vec

 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: /*@C
 58:    VecCheckOrthogonality - Checks (or prints) the level of orthogonality
 59:    of a set of vectors.

 61:    Collective on Vec

 63:    Input parameters:
 64: +  V  - a set of vectors
 65: .  nv - number of V vectors
 66: .  W  - an alternative set of vectors (optional)
 67: .  nw - number of W vectors
 68: .  B  - matrix defining the inner product (optional)
 69: -  viewer - optional visualization context

 71:    Output parameter:
 72: .  lev - level of orthogonality (optional)

 74:    Notes:
 75:    This function computes W'*V and prints the result. It is intended to check
 76:    the level of bi-orthogonality of the vectors in the two sets. If W is equal
 77:    to NULL then V is used, thus checking the orthogonality of the V vectors.

 79:    If matrix B is provided then the check uses the B-inner product, W'*B*V.

 81:    If lev is not NULL, it will contain the maximum entry of matrix
 82:    W'*V - I (in absolute value). Otherwise, the matrix W'*V is printed.

 84:    Level: developer
 85: @*/
 86: PetscErrorCode VecCheckOrthogonality(Vec *V,PetscInt nv,Vec *W,PetscInt nw,Mat B,PetscViewer viewer,PetscReal *lev)
 87: {
 89:   PetscInt       i,j;
 90:   PetscScalar    *vals;
 91:   PetscBool      isascii;
 92:   Vec            w;

 95:   if (nv<=0 || nw<=0) return(0);
 98:   if (!lev) {
 99:     if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)*V));
102:     PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
103:     if (!isascii) return(0);
104:   }

106:   PetscMalloc1(nv,&vals);
107:   if (B) {
108:     VecDuplicate(V[0],&w);
109:   }
110:   if (lev) *lev = 0.0;
111:   for (i=0;i<nw;i++) {
112:     if (B) {
113:       if (W) {
114:         MatMultTranspose(B,W[i],w);
115:       } else {
116:         MatMultTranspose(B,V[i],w);
117:       }
118:     } else {
119:       if (W) w = W[i];
120:       else w = V[i];
121:     }
122:     VecMDot(w,nv,V,vals);
123:     for (j=0;j<nv;j++) {
124:       if (lev) *lev = PetscMax(*lev,PetscAbsScalar((j==i)? (vals[j]-1.0): vals[j]));
125:       else {
126: #if !defined(PETSC_USE_COMPLEX)
127:         PetscViewerASCIIPrintf(viewer," %12g  ",(double)vals[j]);
128: #else
129:         PetscViewerASCIIPrintf(viewer," %12g%+12gi ",(double)PetscRealPart(vals[j]),(double)PetscImaginaryPart(vals[j]));
130: #endif
131:       }
132:     }
133:     if (!lev) { PetscViewerASCIIPrintf(viewer,"\n"); }
134:   }
135:   PetscFree(vals);
136:   if (B) {
137:     VecDestroy(&w);
138:   }
139:   return(0);
140: }

142: /*@
143:    VecDuplicateEmpty - Creates a new vector of the same type as an existing vector,
144:    but without internal array.

146:    Collective on Vec

148:    Input parameters:
149: .  v - a vector to mimic

151:    Output parameter:
152: .  newv - location to put new vector

154:    Note:
155:    This is similar to VecDuplicate(), but the new vector does not have an internal
156:    array, so the intended usage is with VecPlaceArray().

158:    Level: developer
159: @*/
160: PetscErrorCode VecDuplicateEmpty(Vec v,Vec *newv)
161: {
163:   PetscBool      sup,cuda,mpi;
164:   PetscInt       N,nloc,bs;


171:   PetscObjectTypeCompareAny((PetscObject)v,&sup,VECSEQ,VECMPI,VECSEQCUDA,VECMPICUDA,"");
172:   if (!sup) SETERRQ1(PetscObjectComm((PetscObject)v),PETSC_ERR_SUP,"Vector type %s not supported",((PetscObject)v)->type_name);
173:   PetscObjectTypeCompareAny((PetscObject)v,&cuda,VECSEQCUDA,VECMPICUDA,"");
174:   PetscObjectTypeCompareAny((PetscObject)v,&mpi,VECMPI,VECMPICUDA,"");
175:   VecGetLocalSize(v,&nloc);
176:   VecGetSize(v,&N);
177:   VecGetBlockSize(v,&bs);

179:   if (cuda) {
180: #if defined(PETSC_HAVE_VECCUDA)
181:     if (mpi) {
182:       VecCreateMPICUDAWithArray(PetscObjectComm((PetscObject)v),bs,nloc,N,NULL,newv);
183:     } else {
184:       VecCreateSeqCUDAWithArray(PetscObjectComm((PetscObject)v),bs,N,NULL,newv);
185:     }
186: #endif
187:   } else {
188:     if (mpi) {
189:       VecCreateMPIWithArray(PetscObjectComm((PetscObject)v),bs,nloc,N,NULL,newv);
190:     } else {
191:       VecCreateSeqWithArray(PetscObjectComm((PetscObject)v),bs,N,NULL,newv);
192:     }
193:   }
194:   return(0);
195: }