Actual source code: test9.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 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,NULL,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(PetscObjectSetName((PetscObject)G,"G"));
48: PetscCall(MatGetOwnershipRange(G,&Istart,&Iend));
49: for (i=Istart;i<Iend;i++) {
50: col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
51: if (i==0) PetscCall(MatSetValues(G,1,&i,PetscMin(4,n-i),col+1,value+1,INSERT_VALUES));
52: else PetscCall(MatSetValues(G,1,&i,PetscMin(5,n-i+1),col,value,INSERT_VALUES));
53: }
54: PetscCall(MatAssemblyBegin(G,MAT_FINAL_ASSEMBLY));
55: PetscCall(MatAssemblyEnd(G,MAT_FINAL_ASSEMBLY));
56: if (verbose) PetscCall(MatView(G,view));
58: /* Create symmetric matrix B (1-D Laplacian) */
59: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
60: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
61: PetscCall(MatSetFromOptions(B));
62: PetscCall(PetscObjectSetName((PetscObject)B,"B"));
64: PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
65: for (i=Istart;i<Iend;i++) {
66: if (i>0) PetscCall(MatSetValue(B,i,i-1,-1.0,INSERT_VALUES));
67: if (i<n-1) PetscCall(MatSetValue(B,i,i+1,-1.0,INSERT_VALUES));
68: PetscCall(MatSetValue(B,i,i,2.0,INSERT_VALUES));
69: }
70: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
71: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
72: PetscCall(MatCreateVecs(B,&t,NULL));
73: if (verbose) PetscCall(MatView(B,view));
75: /* Create BV object X */
76: PetscCall(BVCreate(PETSC_COMM_WORLD,&X));
77: PetscCall(PetscObjectSetName((PetscObject)X,"X"));
78: PetscCall(BVSetSizesFromVec(X,t,kx+2)); /* two extra columns to test active columns */
79: PetscCall(BVSetFromOptions(X));
81: /* Fill X entries */
82: for (j=0;j<kx+2;j++) {
83: PetscCall(BVGetColumn(X,j,&v));
84: PetscCall(VecSet(v,0.0));
85: for (i=0;i<4;i++) {
86: if (i+j<n) {
87: #if defined(PETSC_USE_COMPLEX)
88: alpha = PetscCMPLX((PetscReal)(3*i+j-2),(PetscReal)(2*i));
89: #else
90: alpha = (PetscReal)(3*i+j-2);
91: #endif
92: PetscCall(VecSetValue(v,i+j,alpha,INSERT_VALUES));
93: }
94: }
95: PetscCall(VecAssemblyBegin(v));
96: PetscCall(VecAssemblyEnd(v));
97: PetscCall(BVRestoreColumn(X,j,&v));
98: }
99: if (verbose) PetscCall(BVView(X,view));
101: /* Duplicate BV object and store Z=G*X */
102: PetscCall(BVDuplicate(X,&Z));
103: PetscCall(PetscObjectSetName((PetscObject)Z,"Z"));
104: PetscCall(BVSetActiveColumns(X,0,kx));
105: PetscCall(BVSetActiveColumns(Z,0,kx));
106: PetscCall(BVMatMult(X,G,Z));
107: PetscCall(BVSetActiveColumns(X,lx,kx));
108: PetscCall(BVSetActiveColumns(Z,lx,kx));
110: /* Create BV object Y */
111: PetscCall(BVCreate(PETSC_COMM_WORLD,&Y));
112: PetscCall(PetscObjectSetName((PetscObject)Y,"Y"));
113: PetscCall(BVSetSizesFromVec(Y,t,ky+1));
114: PetscCall(BVSetFromOptions(Y));
115: PetscCall(BVSetActiveColumns(Y,ly,ky));
117: /* Fill Y entries */
118: for (j=0;j<ky+1;j++) {
119: PetscCall(BVGetColumn(Y,j,&v));
120: #if defined(PETSC_USE_COMPLEX)
121: alpha = PetscCMPLX((PetscReal)(j+1)/4.0,-(PetscReal)j);
122: #else
123: alpha = (PetscReal)(j+1)/4.0;
124: #endif
125: PetscCall(VecSet(v,(PetscScalar)(j+1)/4.0));
126: PetscCall(BVRestoreColumn(Y,j,&v));
127: }
128: if (verbose) PetscCall(BVView(Y,view));
130: /* Test BVMatProject for non-symmetric matrix G */
131: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H0));
132: PetscCall(PetscObjectSetName((PetscObject)H0,"H0"));
133: PetscCall(BVMatProject(X,G,Y,H0));
134: if (verbose) PetscCall(MatView(H0,view));
136: /* Test BVMatProject with previously stored G*X */
137: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H1));
138: PetscCall(PetscObjectSetName((PetscObject)H1,"H1"));
139: PetscCall(BVMatProject(Z,NULL,Y,H1));
140: if (verbose) PetscCall(MatView(H1,view));
142: /* Check that H0 and H1 are equal */
143: PetscCall(MatAXPY(H0,-1.0,H1,SAME_NONZERO_PATTERN));
144: PetscCall(MatNorm(H0,NORM_1,&norm));
145: if (norm<10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||H0-H1|| < 10*eps\n"));
146: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||H0-H1||=%g\n",(double)norm));
147: PetscCall(MatDestroy(&H0));
148: PetscCall(MatDestroy(&H1));
150: /* Test BVMatProject for symmetric matrix B with orthogonal projection */
151: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,kx,kx,NULL,&H0));
152: PetscCall(PetscObjectSetName((PetscObject)H0,"H0"));
153: PetscCall(BVMatProject(X,B,X,H0));
154: if (verbose) PetscCall(MatView(H0,view));
156: /* Repeat previous test with symmetry flag set */
157: PetscCall(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE));
158: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,kx,kx,NULL,&H1));
159: PetscCall(PetscObjectSetName((PetscObject)H1,"H1"));
160: PetscCall(BVMatProject(X,B,X,H1));
161: if (verbose) PetscCall(MatView(H1,view));
163: /* Check that H0 and H1 are equal */
164: PetscCall(MatAXPY(H0,-1.0,H1,SAME_NONZERO_PATTERN));
165: PetscCall(MatNorm(H0,NORM_1,&norm));
166: if (norm<10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||H0-H1|| < 10*eps\n"));
167: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||H0-H1||=%g\n",(double)norm));
168: PetscCall(MatDestroy(&H0));
169: PetscCall(MatDestroy(&H1));
171: PetscCall(BVDestroy(&X));
172: PetscCall(BVDestroy(&Y));
173: PetscCall(BVDestroy(&Z));
174: PetscCall(MatDestroy(&B));
175: PetscCall(MatDestroy(&G));
176: PetscCall(VecDestroy(&t));
177: PetscCall(SlepcFinalize());
178: return 0;
179: }
181: /*TEST
183: testset:
184: output_file: output/test9_1.out
185: test:
186: suffix: 1
187: args: -bv_type {{vecs contiguous svec mat}shared output}
188: test:
189: suffix: 1_svec_vecs
190: args: -bv_type svec -bv_matmult vecs
191: test:
192: suffix: 1_cuda
193: args: -bv_type {{svec mat}} -mat_type aijcusparse
194: requires: cuda
195: test:
196: suffix: 1_hip
197: args: -bv_type {{svec mat}} -mat_type aijhipsparse
198: requires: hip
199: test:
200: suffix: 2
201: nsize: 2
202: args: -bv_type {{vecs contiguous svec mat}shared output}
203: test:
204: suffix: 2_svec_vecs
205: nsize: 2
206: args: -bv_type svec -bv_matmult vecs
208: TEST*/