Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
11 : static char help[] = "Test BV operations with indefinite inner product.\n\n";
12 :
13 : #include <slepcbv.h>
14 :
15 8 : int main(int argc,char **argv)
16 : {
17 8 : Vec t,v,w,omega;
18 8 : Mat B,M;
19 8 : BV X,Y;
20 8 : PetscInt i,j,n=10,k=5,l,Istart,Iend;
21 8 : PetscScalar alpha;
22 8 : PetscReal nrm;
23 8 : PetscViewer view;
24 8 : PetscBool verbose;
25 :
26 8 : PetscFunctionBeginUser;
27 8 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
28 8 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
29 8 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
30 8 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
31 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test BV with indefinite inner product (n=%" PetscInt_FMT ", k=%" PetscInt_FMT ").\n",n,k));
32 :
33 : /* Create inner product matrix (standard involutionary permutation) */
34 8 : PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
35 8 : PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
36 8 : PetscCall(MatSetFromOptions(B));
37 8 : PetscCall(PetscObjectSetName((PetscObject)B,"B"));
38 :
39 8 : PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
40 88 : for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(B,i,n-i-1,1.0,INSERT_VALUES));
41 8 : PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
42 8 : PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
43 8 : PetscCall(MatCreateVecs(B,&t,NULL));
44 :
45 : /* Create BV object X */
46 8 : PetscCall(BVCreate(PETSC_COMM_WORLD,&X));
47 8 : PetscCall(PetscObjectSetName((PetscObject)X,"X"));
48 8 : PetscCall(BVSetSizesFromVec(X,t,k));
49 8 : PetscCall(BVSetFromOptions(X));
50 8 : PetscCall(BVSetMatrix(X,B,PETSC_TRUE));
51 :
52 : /* Set up viewer */
53 8 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view));
54 8 : if (verbose) PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB));
55 :
56 : /* Fill X entries */
57 : l = -3;
58 48 : for (j=0;j<k;j++) {
59 40 : PetscCall(BVGetColumn(X,j,&v));
60 40 : PetscCall(VecSet(v,-1.0));
61 240 : for (i=0;i<n/2;i++) {
62 200 : if (i+j<n) {
63 200 : l = (l + 3*i+j-2) % n;
64 200 : PetscCall(VecSetValue(v,i+j,(PetscScalar)l,INSERT_VALUES));
65 : }
66 : }
67 40 : PetscCall(VecAssemblyBegin(v));
68 40 : PetscCall(VecAssemblyEnd(v));
69 40 : PetscCall(BVRestoreColumn(X,j,&v));
70 : }
71 8 : if (verbose) {
72 0 : PetscCall(MatView(B,view));
73 0 : PetscCall(BVView(X,view));
74 : }
75 :
76 : /* Test BVNormColumn */
77 8 : PetscCall(BVNormColumn(X,0,NORM_2,&nrm));
78 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"B-Norm of X[0] = %g\n",(double)nrm));
79 :
80 : /* Test BVOrthogonalizeColumn */
81 48 : for (j=0;j<k;j++) {
82 40 : PetscCall(BVOrthogonalizeColumn(X,j,NULL,&nrm,NULL));
83 40 : alpha = 1.0/nrm;
84 40 : PetscCall(BVScaleColumn(X,j,alpha));
85 : }
86 8 : if (verbose) PetscCall(BVView(X,view));
87 :
88 : /* Create a copy on Y */
89 8 : PetscCall(BVDuplicate(X,&Y));
90 8 : PetscCall(PetscObjectSetName((PetscObject)Y,"Y"));
91 8 : PetscCall(BVCopy(X,Y));
92 :
93 : /* Check orthogonality */
94 8 : PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M));
95 8 : PetscCall(BVDot(Y,Y,M));
96 8 : PetscCall(VecCreateSeq(PETSC_COMM_SELF,k,&omega));
97 8 : PetscCall(BVGetSignature(Y,omega));
98 8 : PetscCall(VecScale(omega,-1.0));
99 8 : PetscCall(MatDiagonalSet(M,omega,ADD_VALUES));
100 8 : PetscCall(MatNorm(M,NORM_1,&nrm));
101 8 : if (nrm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n"));
102 0 : else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)nrm));
103 :
104 : /* Test BVSetSignature */
105 8 : PetscCall(VecScale(omega,-1.0));
106 8 : PetscCall(BVSetSignature(Y,omega));
107 8 : PetscCall(VecDestroy(&omega));
108 :
109 : /* Test BVApplyMatrix */
110 8 : PetscCall(VecDuplicate(t,&w));
111 8 : PetscCall(BVGetColumn(X,0,&v));
112 8 : PetscCall(BVApplyMatrix(X,v,w));
113 8 : PetscCall(BVApplyMatrix(X,w,t));
114 8 : PetscCall(VecAXPY(t,-1.0,v));
115 8 : PetscCall(BVRestoreColumn(X,0,&v));
116 8 : PetscCall(VecNorm(t,NORM_2,&nrm));
117 8 : PetscCheck(PetscAbsReal(nrm)<10*PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_FP,"Wrong value, nrm = %g",(double)nrm);
118 :
119 8 : PetscCall(BVApplyMatrixBV(X,Y));
120 8 : PetscCall(BVGetColumn(Y,0,&v));
121 8 : PetscCall(VecAXPY(w,-1.0,v));
122 8 : PetscCall(BVRestoreColumn(Y,0,&v));
123 8 : PetscCall(VecNorm(w,NORM_2,&nrm));
124 8 : PetscCheck(PetscAbsReal(nrm)<10*PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_FP,"Wrong value, nrm = %g",(double)nrm);
125 :
126 8 : PetscCall(BVDestroy(&X));
127 8 : PetscCall(BVDestroy(&Y));
128 8 : PetscCall(MatDestroy(&M));
129 8 : PetscCall(MatDestroy(&B));
130 8 : PetscCall(VecDestroy(&w));
131 8 : PetscCall(VecDestroy(&t));
132 8 : PetscCall(SlepcFinalize());
133 : return 0;
134 : }
135 :
136 : /*TEST
137 :
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
163 :
164 : TEST*/
|