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