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