Actual source code: test8.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 orthogonalization with selected columns.\n\n";
13: #include <slepcbv.h>
15: int main(int argc,char **argv)
16: {
18: BV X;
19: Vec v,t,z;
20: PetscInt i,j,n=20,k=8;
21: PetscViewer view;
22: PetscBool verbose,*which;
23: PetscReal norm;
24: PetscScalar alpha,*pz;
26: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
27: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
28: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
29: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
30: PetscPrintf(PETSC_COMM_WORLD,"Test BV orthogonalization with selected columns of length %D.\n",n);
32: /* Create template vector */
33: VecCreate(PETSC_COMM_WORLD,&t);
34: VecSetSizes(t,PETSC_DECIDE,n);
35: VecSetFromOptions(t);
37: /* Create BV object X */
38: BVCreate(PETSC_COMM_WORLD,&X);
39: PetscObjectSetName((PetscObject)X,"X");
40: BVSetSizesFromVec(X,t,k);
41: BVSetOrthogonalization(X,BV_ORTHOG_MGS,BV_ORTHOG_REFINE_IFNEEDED,PETSC_DEFAULT,BV_ORTHOG_BLOCK_GS);
42: BVSetFromOptions(X);
44: /* Set up viewer */
45: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
46: if (verbose) {
47: PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
48: }
50: /* Fill X entries */
51: for (j=0;j<k;j++) {
52: BVGetColumn(X,j,&v);
53: VecSet(v,0.0);
54: for (i=0;i<=n/2;i++) {
55: if (i+j<n) {
56: alpha = (3.0*i+j-2)/(2*(i+j+1));
57: VecSetValue(v,i+j,alpha,INSERT_VALUES);
58: }
59: }
60: VecAssemblyBegin(v);
61: VecAssemblyEnd(v);
62: BVRestoreColumn(X,j,&v);
63: }
64: if (verbose) {
65: BVView(X,view);
66: }
68: /* Orthonormalize first k-1 columns */
69: for (j=0;j<k-1;j++) {
70: BVOrthogonalizeColumn(X,j,NULL,&norm,NULL);
71: alpha = 1.0/norm;
72: BVScaleColumn(X,j,alpha);
73: }
74: if (verbose) {
75: BVView(X,view);
76: }
78: /* Select odd columns and orthogonalize last column against those only */
79: PetscMalloc1(k,&which);
80: for (i=0;i<k;i++) which[i] = (i%2)? PETSC_TRUE: PETSC_FALSE;
81: BVOrthogonalizeSomeColumn(X,k-1,which,NULL,NULL,NULL);
82: PetscFree(which);
83: if (verbose) {
84: BVView(X,view);
85: }
87: /* Check orthogonality */
88: PetscPrintf(PETSC_COMM_WORLD,"Orthogonalization coefficients:\n");
89: VecCreateSeq(PETSC_COMM_SELF,k-1,&z);
90: PetscObjectSetName((PetscObject)z,"z");
91: VecGetArray(z,&pz);
92: BVDotColumn(X,k-1,pz);
93: for (i=0;i<k-1;i++) {
94: if (PetscAbsScalar(pz[i])<5.0*PETSC_MACHINE_EPSILON) pz[i]=0.0;
95: }
96: VecRestoreArray(z,&pz);
97: VecView(z,view);
98: VecDestroy(&z);
100: BVDestroy(&X);
101: VecDestroy(&t);
102: SlepcFinalize();
103: return ierr;
104: }
106: /*TEST
108: test:
109: suffix: 1
110: nsize: 1
111: args: -bv_type {{vecs contiguous svec mat}shared output}
112: output_file: output/test8_1.out
114: test:
115: suffix: 1_cuda
116: nsize: 1
117: args: -bv_type svec -vec_type cuda
118: requires: cuda
119: output_file: output/test8_1.out
121: test:
122: suffix: 2
123: nsize: 1
124: args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_refine never
125: requires: !single
126: output_file: output/test8_1.out
128: test:
129: suffix: 3
130: nsize: 1
131: args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_refine always
132: output_file: output/test8_1.out
134: test:
135: suffix: 4
136: nsize: 1
137: args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type mgs
138: output_file: output/test8_1.out
140: test:
141: suffix: 4_cuda
142: nsize: 1
143: args: -bv_type svec -vec_type cuda -bv_orthog_type mgs
144: requires: cuda
145: output_file: output/test8_1.out
147: TEST*/