Actual source code: test18.c

slepc-3.22.1 2024-10-28
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: static char help[] = "Test BVNormalize().\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 17:   BV             X,Y,Z;
 18:   Mat            B;
 19:   Vec            v,t;
 20:   PetscInt       i,j,n=20,k=8,l=3,Istart,Iend;
 21:   PetscViewer    view;
 22:   PetscBool      verbose;
 23:   PetscReal      norm,error;
 24:   PetscScalar    alpha;
 25: #if !defined(PETSC_USE_COMPLEX)
 26:   PetscScalar    *eigi;
 27:   PetscRandom    rand;
 28:   PetscReal      normr,normi;
 29:   Vec            vi;
 30: #endif

 32:   PetscFunctionBeginUser;
 33:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 34:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 35:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
 36:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL));
 37:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
 38:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test BV normalization with %" PetscInt_FMT " columns of length %" PetscInt_FMT ".\n",k,n));

 40:   /* Create template vector */
 41:   PetscCall(VecCreate(PETSC_COMM_WORLD,&t));
 42:   PetscCall(VecSetSizes(t,PETSC_DECIDE,n));
 43:   PetscCall(VecSetFromOptions(t));

 45:   /* Create BV object X */
 46:   PetscCall(BVCreate(PETSC_COMM_WORLD,&X));
 47:   PetscCall(PetscObjectSetName((PetscObject)X,"X"));
 48:   PetscCall(BVSetSizesFromVec(X,t,k));
 49:   PetscCall(BVSetFromOptions(X));

 51:   /* Set up viewer */
 52:   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view));
 53:   if (verbose) PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB));

 55:   /* Fill X entries */
 56:   for (j=0;j<k;j++) {
 57:     PetscCall(BVGetColumn(X,j,&v));
 58:     PetscCall(VecSet(v,0.0));
 59:     for (i=0;i<=n/2;i++) {
 60:       if (i+j<n) {
 61:         alpha = (3.0*i+j-2)/(2*(i+j+1));
 62:         PetscCall(VecSetValue(v,i+j,alpha,INSERT_VALUES));
 63:       }
 64:     }
 65:     PetscCall(VecAssemblyBegin(v));
 66:     PetscCall(VecAssemblyEnd(v));
 67:     PetscCall(BVRestoreColumn(X,j,&v));
 68:   }
 69:   if (verbose) PetscCall(BVView(X,view));

 71:   /* Create copies on Y and Z */
 72:   PetscCall(BVDuplicate(X,&Y));
 73:   PetscCall(PetscObjectSetName((PetscObject)Y,"Y"));
 74:   PetscCall(BVCopy(X,Y));
 75:   PetscCall(BVDuplicate(X,&Z));
 76:   PetscCall(PetscObjectSetName((PetscObject)Z,"Z"));
 77:   PetscCall(BVCopy(X,Z));
 78:   PetscCall(BVSetActiveColumns(X,l,k));
 79:   PetscCall(BVSetActiveColumns(Y,l,k));
 80:   PetscCall(BVSetActiveColumns(Z,l,k));

 82:   /* Test BVNormalize */
 83:   PetscCall(BVNormalize(X,NULL));
 84:   if (verbose) PetscCall(BVView(X,view));

 86:   /* Check unit norm of columns */
 87:   error = 0.0;
 88:   for (j=l;j<k;j++) {
 89:     PetscCall(BVNormColumn(X,j,NORM_2,&norm));
 90:     error = PetscMax(error,PetscAbsReal(norm-PetscRealConstant(1.0)));
 91:   }
 92:   if (error<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Deviation from normalized vectors < 100*eps\n"));
 93:   else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Deviation from normalized vectors: %g\n",(double)norm));

 95:   /* Create inner product matrix */
 96:   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
 97:   PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
 98:   PetscCall(MatSetFromOptions(B));
 99:   PetscCall(PetscObjectSetName((PetscObject)B,"B"));

101:   PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
102:   for (i=Istart;i<Iend;i++) {
103:     if (i>0) PetscCall(MatSetValue(B,i,i-1,-1.0,INSERT_VALUES));
104:     if (i<n-1) PetscCall(MatSetValue(B,i,i+1,-1.0,INSERT_VALUES));
105:     PetscCall(MatSetValue(B,i,i,2.0,INSERT_VALUES));
106:   }
107:   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
108:   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
109:   if (verbose) PetscCall(MatView(B,view));

111:   /* Test BVNormalize with B-norm */
112:   PetscCall(BVSetMatrix(Y,B,PETSC_FALSE));
113:   PetscCall(BVNormalize(Y,NULL));
114:   if (verbose) PetscCall(BVView(Y,view));

116:   /* Check unit B-norm of columns */
117:   error = 0.0;
118:   for (j=l;j<k;j++) {
119:     PetscCall(BVNormColumn(Y,j,NORM_2,&norm));
120:     error = PetscMax(error,PetscAbsReal(norm-PetscRealConstant(1.0)));
121:   }
122:   if (error<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Deviation from B-normalized vectors < 100*eps\n"));
123:   else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Deviation from B-normalized vectors: %g\n",(double)norm));

125: #if !defined(PETSC_USE_COMPLEX)
126:   /* fill imaginary parts */
127:   PetscCall(PetscCalloc1(k,&eigi));
128:   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
129:   PetscCall(PetscRandomSetFromOptions(rand));
130:   for (j=l+1;j<k-1;j+=5) {
131:     PetscCall(PetscRandomGetValue(rand,&alpha));
132:     eigi[j]   =  alpha;
133:     eigi[j+1] = -alpha;
134:   }
135:   PetscCall(PetscRandomDestroy(&rand));
136:   if (verbose) {
137:     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,k,eigi,&v));
138:     PetscCall(VecView(v,view));
139:     PetscCall(VecDestroy(&v));
140:   }

142:   /* Test BVNormalize with complex conjugate columns */
143:   PetscCall(BVNormalize(Z,eigi));
144:   if (verbose) PetscCall(BVView(Z,view));

146:   /* Check unit norm of (complex conjugate) columns */
147:   error = 0.0;
148:   for (j=l;j<k;j++) {
149:     if (eigi[j]) {
150:       PetscCall(BVGetColumn(Z,j,&v));
151:       PetscCall(BVGetColumn(Z,j+1,&vi));
152:       PetscCall(VecNormBegin(v,NORM_2,&normr));
153:       PetscCall(VecNormBegin(vi,NORM_2,&normi));
154:       PetscCall(VecNormEnd(v,NORM_2,&normr));
155:       PetscCall(VecNormEnd(vi,NORM_2,&normi));
156:       PetscCall(BVRestoreColumn(Z,j+1,&vi));
157:       PetscCall(BVRestoreColumn(Z,j,&v));
158:       norm = SlepcAbsEigenvalue(normr,normi);
159:       j++;
160:     } else {
161:       PetscCall(BVGetColumn(Z,j,&v));
162:       PetscCall(VecNorm(v,NORM_2,&norm));
163:       PetscCall(BVRestoreColumn(Z,j,&v));
164:     }
165:     error = PetscMax(error,PetscAbsReal(norm-PetscRealConstant(1.0)));
166:   }
167:   if (error<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Deviation from normalized conjugate vectors < 100*eps\n"));
168:   else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Deviation from normalized conjugate vectors: %g\n",(double)norm));
169:   PetscCall(PetscFree(eigi));
170: #endif

172:   PetscCall(BVDestroy(&X));
173:   PetscCall(BVDestroy(&Y));
174:   PetscCall(BVDestroy(&Z));
175:   PetscCall(MatDestroy(&B));
176:   PetscCall(VecDestroy(&t));
177:   PetscCall(SlepcFinalize());
178:   return 0;
179: }

181: /*TEST

183:    testset:
184:       args: -n 250 -l 6 -k 15
185:       nsize: {{1 2}}
186:       requires: !complex
187:       output_file: output/test18_1.out
188:       test:
189:          suffix: 1
190:          args: -bv_type {{vecs contiguous svec mat}}
191:       test:
192:          suffix: 1_cuda
193:          args: -bv_type {{svec mat}} -vec_type cuda
194:          requires: cuda
195:       test:
196:          suffix: 1_hip
197:          args: -bv_type {{svec mat}} -vec_type hip
198:          requires: hip

200:    testset:
201:       args: -n 250 -l 6 -k 15
202:       nsize: {{1 2}}
203:       requires: complex
204:       output_file: output/test18_1_complex.out
205:       test:
206:          suffix: 1_complex
207:          args: -bv_type {{vecs contiguous svec mat}}
208:       test:
209:          suffix: 1_cuda_complex
210:          args: -bv_type {{svec mat}} -vec_type cuda
211:          requires: cuda
212:       test:
213:          suffix: 1_hip_complex
214:          args: -bv_type {{svec mat}} -vec_type hip
215:          requires: hip

217: TEST*/