Actual source code: test3.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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: {
18: Vec t,v;
19: Mat B,M;
20: BV X;
21: PetscInt i,j,n=10,k=5,Istart,Iend;
22: PetscScalar alpha;
23: PetscReal nrm;
24: PetscViewer view;
25: PetscBool verbose;
27: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
28: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
29: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
30: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
31: PetscPrintf(PETSC_COMM_WORLD,"Test BV with non-standard inner product (n=%D, k=%D).\n",n,k);
33: /* Create inner product matrix */
34: MatCreate(PETSC_COMM_WORLD,&B);
35: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
36: MatSetFromOptions(B);
37: MatSetUp(B);
38: PetscObjectSetName((PetscObject)B,"B");
40: MatGetOwnershipRange(B,&Istart,&Iend);
41: for (i=Istart;i<Iend;i++) {
42: if (i>0) { MatSetValue(B,i,i-1,-1.0,INSERT_VALUES); }
43: if (i<n-1) { MatSetValue(B,i,i+1,-1.0,INSERT_VALUES); }
44: MatSetValue(B,i,i,2.0,INSERT_VALUES);
45: }
46: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
47: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
48: MatCreateVecs(B,&t,NULL);
50: /* Create BV object X */
51: BVCreate(PETSC_COMM_WORLD,&X);
52: PetscObjectSetName((PetscObject)X,"X");
53: BVSetSizesFromVec(X,t,k);
54: BVSetFromOptions(X);
55: BVSetMatrix(X,B,PETSC_FALSE);
57: /* Set up viewer */
58: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
59: if (verbose) {
60: PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
61: }
63: /* Fill X entries */
64: for (j=0;j<k;j++) {
65: BVGetColumn(X,j,&v);
66: VecSet(v,0.0);
67: for (i=0;i<4;i++) {
68: if (i+j<n) {
69: VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
70: }
71: }
72: VecAssemblyBegin(v);
73: VecAssemblyEnd(v);
74: BVRestoreColumn(X,j,&v);
75: }
76: if (verbose) {
77: MatView(B,view);
78: BVView(X,view);
79: }
81: /* Test BVNormColumn */
82: BVNormColumn(X,0,NORM_2,&nrm);
83: PetscPrintf(PETSC_COMM_WORLD,"B-Norm of X[0] = %g\n",(double)nrm);
85: /* Test BVOrthogonalizeColumn */
86: for (j=0;j<k;j++) {
87: BVOrthogonalizeColumn(X,j,NULL,&nrm,NULL);
88: alpha = 1.0/nrm;
89: BVScaleColumn(X,j,alpha);
90: }
91: if (verbose) {
92: BVView(X,view);
93: }
95: /* Check orthogonality */
96: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
97: BVDot(X,X,M);
98: MatShift(M,-1.0);
99: MatNorm(M,NORM_1,&nrm);
100: if (nrm<100*PETSC_MACHINE_EPSILON) {
101: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
102: } else {
103: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)nrm);
104: }
106: /* Test BVNormVecBegin/End */
107: BVGetColumn(X,0,&v);
108: BVNormVecBegin(X,v,NORM_1,&nrm);
109: BVNormVecEnd(X,v,NORM_1,&nrm);
110: BVRestoreColumn(X,0,&v);
111: PetscPrintf(PETSC_COMM_WORLD,"B-Norm of X[0] = %g\n",(double)nrm);
113: BVDestroy(&X);
114: MatDestroy(&M);
115: MatDestroy(&B);
116: VecDestroy(&t);
117: SlepcFinalize();
118: return ierr;
119: }
121: /*TEST
123: test:
124: suffix: 1
125: nsize: 1
126: args: -bv_type {{vecs contiguous svec mat}shared output}
127: output_file: output/test3_1.out
129: test:
130: suffix: 1_svec_vecs
131: nsize: 1
132: args: -bv_type svec -bv_matmult vecs
133: output_file: output/test3_1.out
135: test:
136: suffix: 1_cuda
137: nsize: 1
138: args: -bv_type svec -mat_type aijcusparse
139: requires: cuda
140: output_file: output/test3_1.out
142: test:
143: suffix: 2
144: nsize: 2
145: args: -bv_type {{vecs contiguous svec mat}shared output}
146: output_file: output/test3_1.out
148: test:
149: suffix: 3
150: nsize: 2
151: args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type mgs
152: output_file: output/test3_1.out
154: TEST*/