Actual source code: test5.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 indefinite inner product.\n\n";
13: #include <slepcbv.h>
15: int main(int argc,char **argv)
16: {
17: Vec t,v,w,omega;
18: Mat B,M;
19: BV X,Y;
20: PetscInt i,j,n=10,k=5,l,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 indefinite inner product (n=%" PetscInt_FMT ", k=%" PetscInt_FMT ").\n",n,k));
33: /* Create inner product matrix (standard involutionary permutation) */
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++) PetscCall(MatSetValue(B,i,n-i-1,1.0,INSERT_VALUES));
41: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
42: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
43: PetscCall(MatCreateVecs(B,&t,NULL));
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));
50: PetscCall(BVSetMatrix(X,B,PETSC_TRUE));
52: /* Set up viewer */
53: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view));
54: if (verbose) PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB));
56: /* Fill X entries */
57: l = -3;
58: for (j=0;j<k;j++) {
59: PetscCall(BVGetColumn(X,j,&v));
60: PetscCall(VecSet(v,-1.0));
61: for (i=0;i<n/2;i++) {
62: if (i+j<n) {
63: l = (l + 3*i+j-2) % n;
64: PetscCall(VecSetValue(v,i+j,(PetscScalar)l,INSERT_VALUES));
65: }
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: /* Create a copy on Y */
89: PetscCall(BVDuplicate(X,&Y));
90: PetscCall(PetscObjectSetName((PetscObject)Y,"Y"));
91: PetscCall(BVCopy(X,Y));
93: /* Check orthogonality */
94: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M));
95: PetscCall(BVDot(Y,Y,M));
96: PetscCall(VecCreateSeq(PETSC_COMM_SELF,k,&omega));
97: PetscCall(BVGetSignature(Y,omega));
98: PetscCall(VecScale(omega,-1.0));
99: PetscCall(MatDiagonalSet(M,omega,ADD_VALUES));
100: PetscCall(MatNorm(M,NORM_1,&nrm));
101: if (nrm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n"));
102: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)nrm));
104: /* Test BVSetSignature */
105: PetscCall(VecScale(omega,-1.0));
106: PetscCall(BVSetSignature(Y,omega));
107: PetscCall(VecDestroy(&omega));
109: /* Test BVApplyMatrix */
110: PetscCall(VecDuplicate(t,&w));
111: PetscCall(BVGetColumn(X,0,&v));
112: PetscCall(BVApplyMatrix(X,v,w));
113: PetscCall(BVApplyMatrix(X,w,t));
114: PetscCall(VecAXPY(t,-1.0,v));
115: PetscCall(BVRestoreColumn(X,0,&v));
116: PetscCall(VecNorm(t,NORM_2,&nrm));
117: PetscCheck(PetscAbsReal(nrm)<10*PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_FP,"Wrong value, nrm = %g",(double)nrm);
119: PetscCall(BVApplyMatrixBV(X,Y));
120: PetscCall(BVGetColumn(Y,0,&v));
121: PetscCall(VecAXPY(w,-1.0,v));
122: PetscCall(BVRestoreColumn(Y,0,&v));
123: PetscCall(VecNorm(w,NORM_2,&nrm));
124: PetscCheck(PetscAbsReal(nrm)<10*PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_FP,"Wrong value, nrm = %g",(double)nrm);
126: PetscCall(BVDestroy(&X));
127: PetscCall(BVDestroy(&Y));
128: PetscCall(MatDestroy(&M));
129: PetscCall(MatDestroy(&B));
130: PetscCall(VecDestroy(&w));
131: PetscCall(VecDestroy(&t));
132: PetscCall(SlepcFinalize());
133: return 0;
134: }
136: /*TEST
138: testset:
139: output_file: output/test5_1.out
140: args: -bv_orthog_refine always
141: test:
142: suffix: 1
143: args: -bv_type {{vecs contiguous svec mat}shared output}
144: test:
145: suffix: 1_cuda
146: args: -bv_type {{svec mat}} -mat_type aijcusparse
147: requires: cuda
148: test:
149: suffix: 1_hip
150: args: -bv_type {{svec mat}} -mat_type aijhipsparse
151: requires: hip
152: test:
153: suffix: 2
154: args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type mgs
155: test:
156: suffix: 2_cuda
157: args: -bv_type {{svec mat}} -mat_type aijcusparse -bv_orthog_type mgs
158: requires: cuda
159: test:
160: suffix: 2_hip
161: args: -bv_type {{svec mat}} -mat_type aijhipsparse -bv_orthog_type mgs
162: requires: hip
164: TEST*/