Actual source code: test20.c
slepc-3.12.1 2019-11-08
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 DSGNHEP with upper quasi-triangular matrix pair.\n\n";
13: #include <slepcds.h>
15: int main(int argc,char **argv)
16: {
18: DS ds;
19: PetscScalar *A,*B,*X;
20: PetscReal rnorm,aux;
21: PetscInt i,j,n=10,ld;
22: PetscViewer viewer;
23: PetscBool verbose;
25: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
26: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
27: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GNHEP - dimension %D.\n",n);
28: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
30: /* Create DS object */
31: DSCreate(PETSC_COMM_WORLD,&ds);
32: DSSetType(ds,DSGNHEP);
33: DSSetFromOptions(ds);
34: ld = n+2; /* test leading dimension larger than n */
35: DSAllocate(ds,ld);
36: DSSetDimensions(ds,n,0,0,0);
38: /* Set up viewer */
39: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
40: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
41: DSView(ds,viewer);
42: PetscViewerPopFormat(viewer);
43: if (verbose) {
44: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
45: }
47: /* Fill A,B with upper quasi-triangular matrices */
48: DSGetArray(ds,DS_MAT_A,&A);
49: DSGetArray(ds,DS_MAT_B,&B);
50: PetscArrayzero(A,ld*n);
51: for (i=0;i<n;i++) A[i+i*ld]=2.0;
52: for (j=1;j<3;j++) {
53: for (i=0;i<n-j;i++) A[i+(i+j)*ld]=0.001;
54: }
55: PetscArrayzero(B,ld*n);
56: for (i=0;i<n;i++) B[i+i*ld]=1.0;
57: B[1+0*ld]=B[0+1*ld]=PETSC_MACHINE_EPSILON;
58: for (i=1;i<n;i+=3) {
59: A[i+(i-1)*ld]=-A[(i-1)+i*ld];
60: }
61: DSRestoreArray(ds,DS_MAT_A,&A);
62: DSRestoreArray(ds,DS_MAT_B,&B);
63: DSSetState(ds,DS_STATE_INTERMEDIATE);
65: if (verbose) {
66: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
67: DSView(ds,viewer);
68: }
70: /* Eigenvectors */
71: j = 0;
72: DSVectors(ds,DS_MAT_X,&j,&rnorm); /* first eigenvector */
73: PetscPrintf(PETSC_COMM_WORLD,"Value of rnorm for 2nd vector = %.3f\n",(double)rnorm);
74: DSVectors(ds,DS_MAT_X,NULL,NULL); /* all eigenvectors */
75: j = 0;
76: rnorm = 0.0;
77: DSGetArray(ds,DS_MAT_X,&X);
78: for (i=0;i<n;i++) {
79: #if defined(PETSC_USE_COMPLEX)
80: aux = PetscAbsScalar(X[i+j*ld]);
81: #else
82: aux = SlepcAbsEigenvalue(X[i+j*ld],X[i+(j+1)*ld]);
83: #endif
84: rnorm += aux*aux;
85: }
86: DSRestoreArray(ds,DS_MAT_X,&X);
87: rnorm = PetscSqrtReal(rnorm);
88: PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st columns = %.3f\n",(double)rnorm);
89: if (verbose) {
90: PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n");
91: DSView(ds,viewer);
92: }
94: DSDestroy(&ds);
95: SlepcFinalize();
96: return ierr;
97: }
99: /*TEST
101: test:
102: suffix: 1
104: TEST*/