Actual source code: test2.c
slepc-3.18.3 2023-03-24
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 orthogonalization functions.\n\n";
13: #include <slepcbv.h>
15: int main(int argc,char **argv)
16: {
17: BV X,Y,Z;
18: Mat M,R;
19: Vec v,t,e;
20: PetscInt i,j,n=20,k=8;
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 orthogonalization with %" PetscInt_FMT " columns of length %" PetscInt_FMT ".\n",k,n);
34: if (condn>1.0) PetscPrintf(PETSC_COMM_WORLD," - Using a random BV 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: VecSetValue(v,i+j,alpha,INSERT_VALUES);
60: }
61: }
62: VecAssemblyBegin(v);
63: VecAssemblyEnd(v);
64: BVRestoreColumn(X,j,&v);
65: }
66: } else BVSetRandomCond(X,condn);
67: if (verbose) BVView(X,view);
69: /* Create copies on Y and Z */
70: BVDuplicate(X,&Y);
71: PetscObjectSetName((PetscObject)Y,"Y");
72: BVCopy(X,Y);
73: BVDuplicate(X,&Z);
74: PetscObjectSetName((PetscObject)Z,"Z");
75: BVCopy(X,Z);
77: /* Test BVOrthogonalizeColumn */
78: for (j=0;j<k;j++) {
79: BVOrthogonalizeColumn(X,j,NULL,&norm,NULL);
80: alpha = 1.0/norm;
81: BVScaleColumn(X,j,alpha);
82: }
83: if (verbose) BVView(X,view);
85: /* Check orthogonality */
86: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
87: BVDot(X,X,M);
88: MatShift(M,-1.0);
89: MatNorm(M,NORM_1,&norm);
90: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
91: else PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);
93: /* Test BVOrthogonalize */
94: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&R);
95: PetscObjectSetName((PetscObject)R,"R");
96: BVOrthogonalize(Y,R);
97: if (verbose) {
98: BVView(Y,view);
99: MatView(R,view);
100: }
102: /* Check orthogonality */
103: BVDot(Y,Y,M);
104: MatShift(M,-1.0);
105: MatNorm(M,NORM_1,&norm);
106: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
107: else PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);
109: /* Check residual */
110: BVMult(Z,-1.0,1.0,Y,R);
111: BVNorm(Z,NORM_FROBENIUS,&norm);
112: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Residual ||X-QR|| < 100*eps\n");
113: else PetscPrintf(PETSC_COMM_WORLD,"Residual ||X-QR||: %g\n",(double)norm);
115: /* Test BVOrthogonalizeVec */
116: VecDuplicate(t,&e);
117: VecSet(e,1.0);
118: BVOrthogonalizeVec(X,e,NULL,&norm,NULL);
119: PetscPrintf(PETSC_COMM_WORLD,"Norm of ones(n,1) after orthogonalizing against X: %g\n",(double)norm);
121: MatDestroy(&M);
122: MatDestroy(&R);
123: BVDestroy(&X);
124: BVDestroy(&Y);
125: BVDestroy(&Z);
126: VecDestroy(&e);
127: VecDestroy(&t);
128: SlepcFinalize();
129: return 0;
130: }
132: /*TEST
134: testset:
135: output_file: output/test2_1.out
136: test:
137: suffix: 1
138: args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type cgs
139: test:
140: suffix: 1_cuda
141: args: -bv_type svec -vec_type cuda -bv_orthog_type cgs
142: requires: cuda
143: test:
144: suffix: 2
145: args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type mgs
146: test:
147: suffix: 2_cuda
148: args: -bv_type svec -vec_type cuda -bv_orthog_type mgs
149: requires: cuda
151: test:
152: suffix: 3
153: nsize: 1
154: args: -bv_type {{vecs contiguous svec mat}shared output} -condn 1e8
155: requires: !single
156: filter: grep -v "against"
157: output_file: output/test2_3.out
159: TEST*/