Actual source code: test17.c
slepc-3.13.0 2020-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: #if defined(PETSC_USE_COMPLEX)
64: alpha += (1.2*i+j-2)/(0.1*(i+j+1))*PETSC_i;
65: #endif
66: VecSetValue(v,i+j,alpha,INSERT_VALUES);
67: }
68: }
69: VecAssemblyBegin(v);
70: VecAssemblyEnd(v);
71: BVRestoreColumn(X,j,&v);
72: }
73: } else {
74: BVSetRandomCond(X,condn);
75: }
76: if (verbose) {
77: BVView(X,view);
78: }
80: /* Create Y and fill its entries */
81: BVDuplicate(X,&Y);
82: PetscObjectSetName((PetscObject)Y,"Y");
83: if (condn==1.0) {
84: for (j=0;j<k;j++) {
85: BVGetColumn(Y,j,&v);
86: VecSet(v,0.0);
87: for (i=PetscMin(n,2+(2*j)%6);i<PetscMin(n,6+(3*j)%9);i++) {
88: if (i%5 && i!=j) {
89: alpha = (1.5*i+j)/(2.2*(i-j));
90: VecSetValue(v,i+j,alpha,INSERT_VALUES);
91: }
92: }
93: VecAssemblyBegin(v);
94: VecAssemblyEnd(v);
95: BVRestoreColumn(Y,j,&v);
96: }
97: } else {
98: BVSetRandomCond(Y,condn);
99: }
100: if (verbose) {
101: BVView(Y,view);
102: }
104: /* Test BVBiorthonormalizeColumn */
105: for (j=0;j<k;j++) {
106: BVBiorthonormalizeColumn(X,Y,j,NULL);
107: }
108: if (verbose) {
109: BVView(X,view);
110: BVView(Y,view);
111: }
113: /* Check orthogonality */
114: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
115: PetscObjectSetName((PetscObject)M,"M");
116: BVDot(X,Y,M);
117: if (verbose) {
118: MatView(M,view);
119: }
120: MatShift(M,-1.0);
121: MatNorm(M,NORM_1,&norm);
122: if (norm<200*PETSC_MACHINE_EPSILON) {
123: PetscPrintf(PETSC_COMM_WORLD,"Level of bi-orthogonality < 200*eps\n");
124: } else {
125: PetscPrintf(PETSC_COMM_WORLD,"Level of bi-orthogonality: %g\n",(double)norm);
126: }
128: MatDestroy(&M);
129: BVDestroy(&X);
130: BVDestroy(&Y);
131: VecDestroy(&t);
132: SlepcFinalize();
133: return ierr;
134: }
136: /*TEST
138: test:
139: suffix: 1
140: nsize: 1
141: args: -bv_type {{vecs contiguous svec mat}} -bv_orthog_type cgs
142: requires: double
143: output_file: output/test17_1.out
145: test:
146: suffix: 1_cuda
147: nsize: 1
148: args: -bv_type svec -vec_type cuda -bv_orthog_type cgs
149: requires: cuda
150: output_file: output/test17_1.out
152: test:
153: suffix: 2
154: nsize: 1
155: args: -bv_type {{vecs contiguous svec mat}} -bv_orthog_type mgs
156: output_file: output/test17_1.out
158: test:
159: suffix: 2_cuda
160: nsize: 1
161: args: -bv_type svec -vec_type cuda -bv_orthog_type mgs
162: requires: cuda
163: output_file: output/test17_1.out
165: TEST*/