Actual source code: test14.c
slepc-3.19.0 2023-03-31
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 created from a dense Mat.\n\n";
13: #include <slepcbv.h>
15: int main(int argc,char **argv)
16: {
17: BV X;
18: Mat A,B,M;
19: PetscInt i,j,n=20,k=8,Istart,Iend;
20: PetscViewer view;
21: PetscBool verbose;
22: PetscReal norm;
23: PetscScalar alpha;
25: PetscFunctionBeginUser;
26: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,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 created from a dense Mat (length %" PetscInt_FMT ", k=%" PetscInt_FMT ").\n",n,k));
32: /* Create dense matrix */
33: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
34: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,k));
35: PetscCall(MatSetType(A,MATDENSE));
36: PetscCall(MatSetUp(A));
37: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
38: for (j=0;j<k;j++) {
39: for (i=0;i<=n/2;i++) {
40: if (i+j<n && i>=Istart && i<Iend) {
41: alpha = (3.0*i+j-2)/(2*(i+j+1));
42: PetscCall(MatSetValue(A,i+j,j,alpha,INSERT_VALUES));
43: }
44: }
45: }
46: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
47: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
49: /* Create BV object X */
50: PetscCall(BVCreateFromMat(A,&X));
51: PetscCall(BVSetFromOptions(X));
52: PetscCall(PetscObjectSetName((PetscObject)X,"X"));
54: /* Set up viewer */
55: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view));
56: if (verbose) {
57: PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB));
58: PetscCall(BVView(X,view));
59: }
61: /* Test BVCreateMat */
62: PetscCall(BVCreateMat(X,&B));
63: PetscCall(MatAXPY(B,-1.0,A,SAME_NONZERO_PATTERN));
64: PetscCall(MatNorm(B,NORM_1,&norm));
65: if (norm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of difference < 100*eps\n"));
66: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of difference: %g\n",(double)norm));
68: /* Test BVOrthogonalize */
69: PetscCall(BVOrthogonalize(X,NULL));
70: if (verbose) PetscCall(BVView(X,view));
72: /* Check orthogonality */
73: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M));
74: PetscCall(MatShift(M,1.0)); /* set leading part to identity */
75: PetscCall(BVDot(X,X,M));
76: PetscCall(MatShift(M,-1.0));
77: PetscCall(MatNorm(M,NORM_1,&norm));
78: if (norm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n"));
79: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm));
81: PetscCall(MatDestroy(&M));
82: PetscCall(MatDestroy(&A));
83: PetscCall(MatDestroy(&B));
84: PetscCall(BVDestroy(&X));
85: PetscCall(SlepcFinalize());
86: return 0;
87: }
89: /*TEST
91: test:
92: suffix: 1
93: nsize: 2
94: args: -bv_type {{vecs contiguous svec mat}shared output}
95: output_file: output/test14_1.out
97: TEST*/