Actual source code: test14.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 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;
26: SlepcInitialize(&argc,&argv,(char*)0,help);
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 created from a dense Mat (length %" PetscInt_FMT ", k=%" PetscInt_FMT ").\n",n,k);
32: /* Create dense matrix */
33: MatCreate(PETSC_COMM_WORLD,&A);
34: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,k);
35: MatSetType(A,MATDENSE);
36: MatSetUp(A);
37: 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: MatSetValue(A,i+j,j,alpha,INSERT_VALUES);
43: }
44: }
45: }
46: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
47: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
49: /* Create BV object X */
50: BVCreateFromMat(A,&X);
51: BVSetFromOptions(X);
52: PetscObjectSetName((PetscObject)X,"X");
54: /* Set up viewer */
55: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
56: if (verbose) {
57: PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
58: BVView(X,view);
59: }
61: /* Test BVCreateMat */
62: BVCreateMat(X,&B);
63: MatAXPY(B,-1.0,A,SAME_NONZERO_PATTERN);
64: MatNorm(B,NORM_1,&norm);
65: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Norm of difference < 100*eps\n");
66: else PetscPrintf(PETSC_COMM_WORLD,"Norm of difference: %g\n",(double)norm);
68: /* Test BVOrthogonalize */
69: BVOrthogonalize(X,NULL);
70: if (verbose) BVView(X,view);
72: /* Check orthogonality */
73: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
74: MatShift(M,1.0); /* set leading part to identity */
75: BVDot(X,X,M);
76: MatShift(M,-1.0);
77: MatNorm(M,NORM_1,&norm);
78: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
79: else PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);
81: MatDestroy(&M);
82: MatDestroy(&A);
83: MatDestroy(&B);
84: BVDestroy(&X);
85: 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*/