Actual source code: test17.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 bi-orthogonalization functions.\n\n";
13: #include <slepcbv.h>
15: int main(int argc,char **argv)
16: {
18: BV X,Y;
19: Mat M;
20: Vec v,t;
21: PetscInt i,j,n=20,k=7;
22: PetscViewer view;
23: PetscBool verbose;
24: PetscReal norm,condn=1.0;
25: PetscScalar alpha;
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: PetscOptionsGetReal(NULL,NULL,"-condn",&condn,NULL);
31: if (condn<1.0) SETERRQ(PETSC_COMM_WORLD,1,"The condition number must be > 1");
32: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
33: PetscPrintf(PETSC_COMM_WORLD,"Test BV bi-orthogonalization with %D columns of length %D.\n",k,n);
34: if (condn>1.0) {
35: PetscPrintf(PETSC_COMM_WORLD," - Using random BVs with condition number = %g\n",(double)condn);
36: }
38: /* Create template vector */
39: VecCreate(PETSC_COMM_WORLD,&t);
40: VecSetSizes(t,PETSC_DECIDE,n);
41: VecSetFromOptions(t);
43: /* Create BV object X */
44: BVCreate(PETSC_COMM_WORLD,&X);
45: PetscObjectSetName((PetscObject)X,"X");
46: BVSetSizesFromVec(X,t,k);
47: BVSetFromOptions(X);
49: /* Set up viewer */
50: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
51: if (verbose) {
52: PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
53: }
55: /* Fill X entries */
56: if (condn==1.0) {
57: for (j=0;j<k;j++) {
58: BVGetColumn(X,j,&v);
59: VecSet(v,0.0);
60: for (i=0;i<=n/2;i++) {
61: if (i+j<n) {
62: alpha = (3.0*i+j-2)/(2*(i+j+1));
63: VecSetValue(v,i+j,alpha,INSERT_VALUES);
64: }
65: }
66: VecAssemblyBegin(v);
67: VecAssemblyEnd(v);
68: BVRestoreColumn(X,j,&v);
69: }
70: } else {
71: BVSetRandomCond(X,condn);
72: }
73: if (verbose) {
74: BVView(X,view);
75: }
77: /* Create Y and fill its entries */
78: BVDuplicate(X,&Y);
79: PetscObjectSetName((PetscObject)Y,"Y");
80: if (condn==1.0) {
81: for (j=0;j<k;j++) {
82: BVGetColumn(Y,j,&v);
83: VecSet(v,0.0);
84: for (i=PetscMin(n,2+(2*j)%6);i<PetscMin(n,6+(3*j)%9);i++) {
85: if (i%5 && i!=j) {
86: alpha = (1.5*i+j)/(2.2*(i-j));
87: VecSetValue(v,i+j,alpha,INSERT_VALUES);
88: }
89: }
90: VecAssemblyBegin(v);
91: VecAssemblyEnd(v);
92: BVRestoreColumn(Y,j,&v);
93: }
94: } else {
95: BVSetRandomCond(Y,condn);
96: }
97: if (verbose) {
98: BVView(Y,view);
99: }
101: /* Test BVBiorthonormalizeColumn */
102: for (j=0;j<k;j++) {
103: BVBiorthonormalizeColumn(X,Y,j,NULL);
104: }
105: if (verbose) {
106: BVView(X,view);
107: BVView(Y,view);
108: }
110: /* Check orthogonality */
111: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
112: PetscObjectSetName((PetscObject)M,"M");
113: BVDot(X,Y,M);
114: if (verbose) {
115: MatView(M,view);
116: }
117: MatShift(M,-1.0);
118: MatNorm(M,NORM_1,&norm);
119: if (norm<100*PETSC_MACHINE_EPSILON) {
120: PetscPrintf(PETSC_COMM_WORLD,"Level of bi-orthogonality < 100*eps\n");
121: } else {
122: PetscPrintf(PETSC_COMM_WORLD,"Level of bi-orthogonality: %g\n",(double)norm);
123: }
125: MatDestroy(&M);
126: BVDestroy(&X);
127: BVDestroy(&Y);
128: VecDestroy(&t);
129: SlepcFinalize();
130: return ierr;
131: }
133: /*TEST
135: test:
136: suffix: 1
137: nsize: 1
138: args: -bv_type {{vecs contiguous svec mat}} -bv_orthog_type cgs
139: requires: double
140: output_file: output/test17_1.out
142: test:
143: suffix: 1_cuda
144: nsize: 1
145: args: -bv_type svec -vec_type cuda -bv_orthog_type cgs
146: requires: cuda
147: output_file: output/test17_1.out
149: test:
150: suffix: 2
151: nsize: 1
152: args: -bv_type {{vecs contiguous svec mat}} -bv_orthog_type mgs
153: output_file: output/test17_1.out
155: test:
156: suffix: 2_cuda
157: nsize: 1
158: args: -bv_type svec -vec_type cuda -bv_orthog_type mgs
159: requires: cuda
160: output_file: output/test17_1.out
162: TEST*/