Actual source code: test17.c
slepc-3.18.2 2023-01-26
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 bi-orthogonalization functions.\n\n";
13: #include <slepcbv.h>
15: int main(int argc,char **argv)
16: {
17: BV X,Y;
18: Mat M;
19: Vec v,t;
20: PetscInt i,j,n=20,k=7;
21: PetscViewer view;
22: PetscBool verbose;
23: PetscReal norm,condn=1.0;
24: PetscScalar alpha;
27: SlepcInitialize(&argc,&argv,(char*)0,help);
28: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
29: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
30: PetscOptionsGetReal(NULL,NULL,"-condn",&condn,NULL);
32: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
33: PetscPrintf(PETSC_COMM_WORLD,"Test BV bi-orthogonalization with %" PetscInt_FMT " columns of length %" PetscInt_FMT ".\n",k,n);
34: if (condn>1.0) PetscPrintf(PETSC_COMM_WORLD," - Using random BVs with condition number = %g\n",(double)condn);
36: /* Create template vector */
37: VecCreate(PETSC_COMM_WORLD,&t);
38: VecSetSizes(t,PETSC_DECIDE,n);
39: VecSetFromOptions(t);
41: /* Create BV object X */
42: BVCreate(PETSC_COMM_WORLD,&X);
43: PetscObjectSetName((PetscObject)X,"X");
44: BVSetSizesFromVec(X,t,k);
45: BVSetFromOptions(X);
47: /* Set up viewer */
48: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
49: if (verbose) PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
51: /* Fill X entries */
52: if (condn==1.0) {
53: for (j=0;j<k;j++) {
54: BVGetColumn(X,j,&v);
55: VecSet(v,0.0);
56: for (i=0;i<=n/2;i++) {
57: if (i+j<n) {
58: alpha = (3.0*i+j-2)/(2*(i+j+1));
59: #if defined(PETSC_USE_COMPLEX)
60: alpha += (1.2*i+j-2)/(0.1*(i+j+1))*PETSC_i;
61: #endif
62: VecSetValue(v,i+j,alpha,INSERT_VALUES);
63: }
64: }
65: VecAssemblyBegin(v);
66: VecAssemblyEnd(v);
67: BVRestoreColumn(X,j,&v);
68: }
69: } else BVSetRandomCond(X,condn);
70: if (verbose) BVView(X,view);
72: /* Create Y and fill its entries */
73: BVDuplicate(X,&Y);
74: PetscObjectSetName((PetscObject)Y,"Y");
75: if (condn==1.0) {
76: for (j=0;j<k;j++) {
77: BVGetColumn(Y,j,&v);
78: VecSet(v,0.0);
79: for (i=PetscMin(n,2+(2*j)%6);i<PetscMin(n,6+(3*j)%9);i++) {
80: if (i%5 && i!=j) {
81: alpha = (1.5*i+j)/(2.2*(i-j));
82: VecSetValue(v,i+j,alpha,INSERT_VALUES);
83: }
84: }
85: VecAssemblyBegin(v);
86: VecAssemblyEnd(v);
87: BVRestoreColumn(Y,j,&v);
88: }
89: } else BVSetRandomCond(Y,condn);
90: if (verbose) BVView(Y,view);
92: /* Test BVBiorthonormalizeColumn */
93: for (j=0;j<k;j++) BVBiorthonormalizeColumn(X,Y,j,NULL);
94: if (verbose) {
95: BVView(X,view);
96: BVView(Y,view);
97: }
99: /* Check orthogonality */
100: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
101: PetscObjectSetName((PetscObject)M,"M");
102: BVDot(X,Y,M);
103: if (verbose) MatView(M,view);
104: MatShift(M,-1.0);
105: MatNorm(M,NORM_1,&norm);
106: if (norm<200*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Level of bi-orthogonality < 200*eps\n");
107: else PetscPrintf(PETSC_COMM_WORLD,"Level of bi-orthogonality: %g\n",(double)norm);
109: MatDestroy(&M);
110: BVDestroy(&X);
111: BVDestroy(&Y);
112: VecDestroy(&t);
113: SlepcFinalize();
114: return 0;
115: }
117: /*TEST
119: testset:
120: output_file: output/test17_1.out
121: test:
122: suffix: 1
123: args: -bv_type {{vecs contiguous svec mat}} -bv_orthog_type cgs
124: requires: double
125: test:
126: suffix: 1_cuda
127: args: -bv_type svec -vec_type cuda -bv_orthog_type cgs
128: requires: cuda
129: test:
130: suffix: 2
131: args: -bv_type {{vecs contiguous svec mat}} -bv_orthog_type mgs
132: test:
133: suffix: 2_cuda
134: args: -bv_type svec -vec_type cuda -bv_orthog_type mgs
135: requires: cuda
137: TEST*/