Actual source code: test9.c
slepc-3.19.2 2023-09-05
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 matrix projection.\n\n";
13: #include <slepcbv.h>
15: int main(int argc,char **argv)
16: {
17: Vec t,v;
18: Mat B,G,H0,H1;
19: BV X,Y,Z;
20: PetscInt i,j,n=20,kx=6,lx=3,ky=5,ly=2,Istart,Iend,col[5];
21: PetscScalar alpha,value[] = { -1, 1, 1, 1, 1 };
22: PetscViewer view;
23: PetscReal norm;
24: PetscBool verbose;
26: PetscFunctionBeginUser;
27: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
28: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
29: PetscCall(PetscOptionsGetInt(NULL,NULL,"-kx",&kx,NULL));
30: PetscCall(PetscOptionsGetInt(NULL,NULL,"-lx",&lx,NULL));
31: PetscCall(PetscOptionsGetInt(NULL,NULL,"-ky",&ky,NULL));
32: PetscCall(PetscOptionsGetInt(NULL,NULL,"-ly",&ly,NULL));
33: PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
34: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test BV projection (n=%" PetscInt_FMT ").\n",n));
35: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"X has %" PetscInt_FMT " active columns (%" PetscInt_FMT " leading columns).\n",kx,lx));
36: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Y has %" PetscInt_FMT " active columns (%" PetscInt_FMT " leading columns).\n",ky,ly));
38: /* Set up viewer */
39: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view));
40: if (verbose) PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB));
42: /* Create non-symmetric matrix G (Toeplitz) */
43: PetscCall(MatCreate(PETSC_COMM_WORLD,&G));
44: PetscCall(MatSetSizes(G,PETSC_DECIDE,PETSC_DECIDE,n,n));
45: PetscCall(MatSetFromOptions(G));
46: PetscCall(MatSetUp(G));
47: PetscCall(PetscObjectSetName((PetscObject)G,"G"));
49: PetscCall(MatGetOwnershipRange(G,&Istart,&Iend));
50: for (i=Istart;i<Iend;i++) {
51: col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
52: if (i==0) PetscCall(MatSetValues(G,1,&i,PetscMin(4,n-i),col+1,value+1,INSERT_VALUES));
53: else PetscCall(MatSetValues(G,1,&i,PetscMin(5,n-i+1),col,value,INSERT_VALUES));
54: }
55: PetscCall(MatAssemblyBegin(G,MAT_FINAL_ASSEMBLY));
56: PetscCall(MatAssemblyEnd(G,MAT_FINAL_ASSEMBLY));
57: if (verbose) PetscCall(MatView(G,view));
59: /* Create symmetric matrix B (1-D Laplacian) */
60: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
61: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
62: PetscCall(MatSetFromOptions(B));
63: PetscCall(MatSetUp(B));
64: PetscCall(PetscObjectSetName((PetscObject)B,"B"));
66: PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
67: for (i=Istart;i<Iend;i++) {
68: if (i>0) PetscCall(MatSetValue(B,i,i-1,-1.0,INSERT_VALUES));
69: if (i<n-1) PetscCall(MatSetValue(B,i,i+1,-1.0,INSERT_VALUES));
70: PetscCall(MatSetValue(B,i,i,2.0,INSERT_VALUES));
71: }
72: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
73: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
74: PetscCall(MatCreateVecs(B,&t,NULL));
75: if (verbose) PetscCall(MatView(B,view));
77: /* Create BV object X */
78: PetscCall(BVCreate(PETSC_COMM_WORLD,&X));
79: PetscCall(PetscObjectSetName((PetscObject)X,"X"));
80: PetscCall(BVSetSizesFromVec(X,t,kx+2)); /* two extra columns to test active columns */
81: PetscCall(BVSetFromOptions(X));
83: /* Fill X entries */
84: for (j=0;j<kx+2;j++) {
85: PetscCall(BVGetColumn(X,j,&v));
86: PetscCall(VecSet(v,0.0));
87: for (i=0;i<4;i++) {
88: if (i+j<n) {
89: #if defined(PETSC_USE_COMPLEX)
90: alpha = PetscCMPLX((PetscReal)(3*i+j-2),(PetscReal)(2*i));
91: #else
92: alpha = (PetscReal)(3*i+j-2);
93: #endif
94: PetscCall(VecSetValue(v,i+j,alpha,INSERT_VALUES));
95: }
96: }
97: PetscCall(VecAssemblyBegin(v));
98: PetscCall(VecAssemblyEnd(v));
99: PetscCall(BVRestoreColumn(X,j,&v));
100: }
101: if (verbose) PetscCall(BVView(X,view));
103: /* Duplicate BV object and store Z=G*X */
104: PetscCall(BVDuplicate(X,&Z));
105: PetscCall(PetscObjectSetName((PetscObject)Z,"Z"));
106: PetscCall(BVSetActiveColumns(X,0,kx));
107: PetscCall(BVSetActiveColumns(Z,0,kx));
108: PetscCall(BVMatMult(X,G,Z));
109: PetscCall(BVSetActiveColumns(X,lx,kx));
110: PetscCall(BVSetActiveColumns(Z,lx,kx));
112: /* Create BV object Y */
113: PetscCall(BVCreate(PETSC_COMM_WORLD,&Y));
114: PetscCall(PetscObjectSetName((PetscObject)Y,"Y"));
115: PetscCall(BVSetSizesFromVec(Y,t,ky+1));
116: PetscCall(BVSetFromOptions(Y));
117: PetscCall(BVSetActiveColumns(Y,ly,ky));
119: /* Fill Y entries */
120: for (j=0;j<ky+1;j++) {
121: PetscCall(BVGetColumn(Y,j,&v));
122: #if defined(PETSC_USE_COMPLEX)
123: alpha = PetscCMPLX((PetscReal)(j+1)/4.0,-(PetscReal)j);
124: #else
125: alpha = (PetscReal)(j+1)/4.0;
126: #endif
127: PetscCall(VecSet(v,(PetscScalar)(j+1)/4.0));
128: PetscCall(BVRestoreColumn(Y,j,&v));
129: }
130: if (verbose) PetscCall(BVView(Y,view));
132: /* Test BVMatProject for non-symmetric matrix G */
133: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H0));
134: PetscCall(PetscObjectSetName((PetscObject)H0,"H0"));
135: PetscCall(BVMatProject(X,G,Y,H0));
136: if (verbose) PetscCall(MatView(H0,view));
138: /* Test BVMatProject with previously stored G*X */
139: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H1));
140: PetscCall(PetscObjectSetName((PetscObject)H1,"H1"));
141: PetscCall(BVMatProject(Z,NULL,Y,H1));
142: if (verbose) PetscCall(MatView(H1,view));
144: /* Check that H0 and H1 are equal */
145: PetscCall(MatAXPY(H0,-1.0,H1,SAME_NONZERO_PATTERN));
146: PetscCall(MatNorm(H0,NORM_1,&norm));
147: if (norm<10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||H0-H1|| < 10*eps\n"));
148: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||H0-H1||=%g\n",(double)norm));
149: PetscCall(MatDestroy(&H0));
150: PetscCall(MatDestroy(&H1));
152: /* Test BVMatProject for symmetric matrix B with orthogonal projection */
153: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,kx,kx,NULL,&H0));
154: PetscCall(PetscObjectSetName((PetscObject)H0,"H0"));
155: PetscCall(BVMatProject(X,B,X,H0));
156: if (verbose) PetscCall(MatView(H0,view));
158: /* Repeat previous test with symmetry flag set */
159: PetscCall(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE));
160: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,kx,kx,NULL,&H1));
161: PetscCall(PetscObjectSetName((PetscObject)H1,"H1"));
162: PetscCall(BVMatProject(X,B,X,H1));
163: if (verbose) PetscCall(MatView(H1,view));
165: /* Check that H0 and H1 are equal */
166: PetscCall(MatAXPY(H0,-1.0,H1,SAME_NONZERO_PATTERN));
167: PetscCall(MatNorm(H0,NORM_1,&norm));
168: if (norm<10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||H0-H1|| < 10*eps\n"));
169: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||H0-H1||=%g\n",(double)norm));
170: PetscCall(MatDestroy(&H0));
171: PetscCall(MatDestroy(&H1));
173: PetscCall(BVDestroy(&X));
174: PetscCall(BVDestroy(&Y));
175: PetscCall(BVDestroy(&Z));
176: PetscCall(MatDestroy(&B));
177: PetscCall(MatDestroy(&G));
178: PetscCall(VecDestroy(&t));
179: PetscCall(SlepcFinalize());
180: return 0;
181: }
183: /*TEST
185: testset:
186: output_file: output/test9_1.out
187: test:
188: suffix: 1
189: args: -bv_type {{vecs contiguous svec mat}shared output}
190: test:
191: suffix: 1_svec_vecs
192: args: -bv_type svec -bv_matmult vecs
193: test:
194: suffix: 1_cuda
195: args: -bv_type svec -mat_type aijcusparse
196: requires: cuda
197: test:
198: suffix: 2
199: nsize: 2
200: args: -bv_type {{vecs contiguous svec mat}shared output}
201: test:
202: suffix: 2_svec_vecs
203: nsize: 2
204: args: -bv_type svec -bv_matmult vecs
206: TEST*/