Actual source code: test4.c
slepc-3.17.0 2022-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 DSGNHEP.\n\n";
13: #include <slepcds.h>
15: int main(int argc,char **argv)
16: {
17: DS ds;
18: SlepcSC sc;
19: PetscScalar *A,*B,*X,*wr,*wi;
20: PetscReal re,im,rnorm,aux;
21: PetscInt i,j,n=10,ld;
22: PetscViewer viewer;
23: PetscBool verbose;
25: SlepcInitialize(&argc,&argv,(char*)0,help);
26: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
27: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GNHEP - dimension %" PetscInt_FMT ".\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);
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) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
45: /* Fill A with Grcar matrix */
46: DSGetArray(ds,DS_MAT_A,&A);
47: PetscArrayzero(A,ld*n);
48: for (i=1;i<n;i++) A[i+(i-1)*ld]=-1.0;
49: for (j=0;j<4;j++) {
50: for (i=0;i<n-j;i++) A[i+(i+j)*ld]=1.0;
51: }
52: DSRestoreArray(ds,DS_MAT_A,&A);
53: /* Fill B with an upper triangular matrix */
54: DSGetArray(ds,DS_MAT_B,&B);
55: PetscArrayzero(B,ld*n);
56: B[0+0*ld]=-1.0;
57: B[0+1*ld]=2.0;
58: for (i=1;i<n;i++) B[i+i*ld]=1.0;
59: DSRestoreArray(ds,DS_MAT_B,&B);
61: if (verbose) {
62: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
63: DSView(ds,viewer);
64: }
66: /* Solve */
67: PetscMalloc2(n,&wr,n,&wi);
68: DSGetSlepcSC(ds,&sc);
69: sc->comparison = SlepcCompareLargestMagnitude;
70: sc->comparisonctx = NULL;
71: sc->map = NULL;
72: sc->mapobj = NULL;
73: DSSolve(ds,wr,wi);
74: DSSort(ds,wr,wi,NULL,NULL,NULL);
75: if (verbose) {
76: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
77: DSView(ds,viewer);
78: }
80: /* Print eigenvalues */
81: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
82: for (i=0;i<n;i++) {
83: #if defined(PETSC_USE_COMPLEX)
84: re = PetscRealPart(wr[i]);
85: im = PetscImaginaryPart(wr[i]);
86: #else
87: re = wr[i];
88: im = wi[i];
89: #endif
90: if (PetscAbs(im)<1e-10) PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re);
91: else PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im);
92: }
94: /* Eigenvectors */
95: j = 1;
96: DSVectors(ds,DS_MAT_X,&j,&rnorm); /* second eigenvector */
97: PetscPrintf(PETSC_COMM_WORLD,"Value of rnorm for 2nd vector = %.3f\n",(double)rnorm);
98: DSVectors(ds,DS_MAT_X,NULL,NULL); /* all eigenvectors */
99: j = 0;
100: rnorm = 0.0;
101: DSGetArray(ds,DS_MAT_X,&X);
102: for (i=0;i<n;i++) {
103: #if defined(PETSC_USE_COMPLEX)
104: aux = PetscAbsScalar(X[i+j*ld]);
105: #else
106: if (PetscAbs(wi[j])==0.0) aux = PetscAbsScalar(X[i+j*ld]);
107: else aux = SlepcAbsEigenvalue(X[i+j*ld],X[i+(j+1)*ld]);
108: #endif
109: rnorm += aux*aux;
110: }
111: DSRestoreArray(ds,DS_MAT_X,&X);
112: rnorm = PetscSqrtReal(rnorm);
113: PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st vector = %.3f\n",(double)rnorm);
114: if (verbose) {
115: PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n");
116: DSView(ds,viewer);
117: }
119: PetscFree2(wr,wi);
120: DSDestroy(&ds);
121: SlepcFinalize();
122: return 0;
123: }
125: /*TEST
127: test:
128: suffix: 1
129: filter: sed -e "s/[+-]\([0-9]\.[0-9]*i\)/+-\\1/" | sed -e "s/+-0\.0*i//"
131: TEST*/