Actual source code: test3.c
slepc-3.22.1 2024-10-28
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 BV operations with non-standard inner product.\n\n";
13: #include <slepcbv.h>
15: int main(int argc,char **argv)
16: {
17: Vec t,v;
18: Mat B,M;
19: BV X;
20: PetscInt i,j,n=10,k=5,Istart,Iend;
21: PetscScalar alpha;
22: PetscReal nrm;
23: PetscViewer view;
24: PetscBool verbose;
26: PetscFunctionBeginUser;
27: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
28: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
29: PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
30: PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
31: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test BV with non-standard inner product (n=%" PetscInt_FMT ", k=%" PetscInt_FMT ").\n",n,k));
33: /* Create inner product matrix */
34: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
35: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
36: PetscCall(MatSetFromOptions(B));
37: PetscCall(PetscObjectSetName((PetscObject)B,"B"));
39: PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
40: for (i=Istart;i<Iend;i++) {
41: if (i>0) PetscCall(MatSetValue(B,i,i-1,-1.0,INSERT_VALUES));
42: if (i<n-1) PetscCall(MatSetValue(B,i,i+1,-1.0,INSERT_VALUES));
43: PetscCall(MatSetValue(B,i,i,2.0,INSERT_VALUES));
44: }
45: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
46: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
47: PetscCall(MatCreateVecs(B,&t,NULL));
49: /* Create BV object X */
50: PetscCall(BVCreate(PETSC_COMM_WORLD,&X));
51: PetscCall(PetscObjectSetName((PetscObject)X,"X"));
52: PetscCall(BVSetSizesFromVec(X,t,k));
53: PetscCall(BVSetFromOptions(X));
54: PetscCall(BVSetMatrix(X,B,PETSC_FALSE));
56: /* Set up viewer */
57: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view));
58: if (verbose) PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB));
60: /* Fill X entries */
61: for (j=0;j<k;j++) {
62: PetscCall(BVGetColumn(X,j,&v));
63: PetscCall(VecSet(v,0.0));
64: for (i=0;i<4;i++) {
65: if (i+j<n) PetscCall(VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES));
66: }
67: PetscCall(VecAssemblyBegin(v));
68: PetscCall(VecAssemblyEnd(v));
69: PetscCall(BVRestoreColumn(X,j,&v));
70: }
71: if (verbose) {
72: PetscCall(MatView(B,view));
73: PetscCall(BVView(X,view));
74: }
76: /* Test BVNormColumn */
77: PetscCall(BVNormColumn(X,0,NORM_2,&nrm));
78: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"B-Norm of X[0] = %g\n",(double)nrm));
80: /* Test BVOrthogonalizeColumn */
81: for (j=0;j<k;j++) {
82: PetscCall(BVOrthogonalizeColumn(X,j,NULL,&nrm,NULL));
83: alpha = 1.0/nrm;
84: PetscCall(BVScaleColumn(X,j,alpha));
85: }
86: if (verbose) PetscCall(BVView(X,view));
88: /* Check orthogonality */
89: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M));
90: PetscCall(BVDot(X,X,M));
91: PetscCall(MatShift(M,-1.0));
92: PetscCall(MatNorm(M,NORM_1,&nrm));
93: if (nrm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n"));
94: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)nrm));
96: /* Test BVNormVecBegin/End */
97: PetscCall(BVGetColumn(X,0,&v));
98: PetscCall(BVNormVecBegin(X,v,NORM_1,&nrm));
99: PetscCall(BVNormVecEnd(X,v,NORM_1,&nrm));
100: PetscCall(BVRestoreColumn(X,0,&v));
101: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"B-Norm of X[0] = %g\n",(double)nrm));
103: PetscCall(BVDestroy(&X));
104: PetscCall(MatDestroy(&M));
105: PetscCall(MatDestroy(&B));
106: PetscCall(VecDestroy(&t));
107: PetscCall(SlepcFinalize());
108: return 0;
109: }
111: /*TEST
113: testset:
114: output_file: output/test3_1.out
115: test:
116: suffix: 1
117: args: -bv_type {{vecs contiguous svec mat}shared output}
118: test:
119: suffix: 1_svec_vecs
120: args: -bv_type svec -bv_matmult vecs
121: test:
122: suffix: 1_cuda
123: args: -bv_type {{svec mat}} -mat_type aijcusparse
124: requires: cuda
125: test:
126: suffix: 1_hip
127: args: -bv_type {{svec mat}} -mat_type aijhipsparse
128: requires: hip
129: test:
130: suffix: 2
131: nsize: 2
132: args: -bv_type {{vecs contiguous svec mat}shared output}
133: test:
134: suffix: 3
135: nsize: 2
136: args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type mgs
138: TEST*/