Actual source code: test1.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 DSNHEP.\n\n";
13: #include <slepcds.h>
15: int main(int argc,char **argv)
16: {
17: DS ds;
18: SlepcSC sc;
19: DSType type;
20: DSStateType state;
21: PetscScalar *A,*X,*Q,*wr,*wi,d;
22: PetscReal re,im,rnorm,aux;
23: PetscInt i,j,n=10,ld,method;
24: PetscViewer viewer;
25: PetscBool verbose,extrarow;
27: SlepcInitialize(&argc,&argv,(char*)0,help);
28: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
29: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NHEP - dimension %" PetscInt_FMT ".\n",n);
30: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
31: PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow);
33: /* Create DS object */
34: DSCreate(PETSC_COMM_WORLD,&ds);
35: DSSetType(ds,DSNHEP);
36: DSSetFromOptions(ds);
37: ld = n+2; /* test leading dimension larger than n */
38: DSAllocate(ds,ld);
39: DSSetDimensions(ds,n,0,0);
40: DSSetExtraRow(ds,extrarow);
42: /* Set up viewer */
43: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
44: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
45: DSView(ds,viewer);
46: PetscViewerPopFormat(viewer);
47: if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
49: /* Fill with Grcar matrix */
50: DSGetArray(ds,DS_MAT_A,&A);
51: for (i=1;i<n;i++) A[i+(i-1)*ld]=-1.0;
52: for (j=0;j<4;j++) {
53: for (i=0;i<n-j;i++) A[i+(i+j)*ld]=1.0;
54: }
55: if (extrarow) A[n+(n-1)*ld]=-1.0;
56: DSRestoreArray(ds,DS_MAT_A,&A);
57: DSSetState(ds,DS_STATE_INTERMEDIATE);
58: if (verbose) {
59: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
60: DSView(ds,viewer);
61: }
63: /* Solve */
64: PetscMalloc2(n,&wr,n,&wi);
65: DSGetSlepcSC(ds,&sc);
66: sc->comparison = SlepcCompareLargestMagnitude;
67: sc->comparisonctx = NULL;
68: sc->map = NULL;
69: sc->mapobj = NULL;
70: DSSolve(ds,wr,wi);
71: DSSort(ds,wr,wi,NULL,NULL,NULL);
72: if (extrarow) DSUpdateExtraRow(ds);
74: DSGetType(ds,&type);
75: DSGetMethod(ds,&method);
76: PetscPrintf(PETSC_COMM_WORLD,"DS of type %s, method used=%" PetscInt_FMT "\n",type,method);
77: DSGetState(ds,&state);
78: PetscPrintf(PETSC_COMM_WORLD,"State after solve: %s\n",DSStateTypes[state]);
80: if (verbose) {
81: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
82: DSView(ds,viewer);
83: }
85: /* Print eigenvalues */
86: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
87: for (i=0;i<n;i++) {
88: #if defined(PETSC_USE_COMPLEX)
89: re = PetscRealPart(wr[i]);
90: im = PetscImaginaryPart(wr[i]);
91: #else
92: re = wr[i];
93: im = wi[i];
94: #endif
95: if (PetscAbs(im)<1e-10) PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re);
96: else PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im);
97: }
99: if (extrarow) {
100: /* Check that extra row is correct */
101: DSGetArray(ds,DS_MAT_A,&A);
102: DSGetArray(ds,DS_MAT_Q,&Q);
103: d = 0.0;
104: for (i=0;i<n;i++) d += A[n+i*ld]+Q[n-1+i*ld];
105: if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in the extra row of %g\n",(double)PetscAbsScalar(d));
106: DSRestoreArray(ds,DS_MAT_A,&A);
107: DSRestoreArray(ds,DS_MAT_Q,&Q);
108: }
110: /* Eigenvectors */
111: j = 2;
112: DSVectors(ds,DS_MAT_X,&j,&rnorm); /* third eigenvector */
113: PetscPrintf(PETSC_COMM_WORLD,"Value of rnorm for 3rd vector = %.3f\n",(double)rnorm);
114: DSVectors(ds,DS_MAT_X,NULL,NULL); /* all eigenvectors */
115: j = 0;
116: rnorm = 0.0;
117: DSGetArray(ds,DS_MAT_X,&X);
118: for (i=0;i<n;i++) {
119: #if defined(PETSC_USE_COMPLEX)
120: aux = PetscAbsScalar(X[i+j*ld]);
121: #else
122: if (PetscAbs(wi[j])==0.0) aux = PetscAbsScalar(X[i+j*ld]);
123: else aux = SlepcAbsEigenvalue(X[i+j*ld],X[i+(j+1)*ld]);
124: #endif
125: rnorm += aux*aux;
126: }
127: DSRestoreArray(ds,DS_MAT_X,&X);
128: rnorm = PetscSqrtReal(rnorm);
129: PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st vector = %.3f\n",(double)rnorm);
130: if (verbose) {
131: PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n");
132: DSView(ds,viewer);
133: }
135: PetscFree2(wr,wi);
136: DSDestroy(&ds);
137: SlepcFinalize();
138: return 0;
139: }
141: /*TEST
143: testset:
144: filter: sed -e "s/[+-]\([0-9]\.[0-9]*i\)/+-\\1/" | sed -e "s/extrarow//"
145: output_file: output/test1_1.out
146: requires: !single
147: test:
148: suffix: 1
149: test:
150: suffix: 2
151: args: -extrarow
153: TEST*/