Actual source code: test6.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 DSGHIEP with compact storage.\n\n";
13: #include <slepcds.h>
15: int main(int argc,char **argv)
16: {
17: DS ds;
18: SlepcSC sc;
19: PetscReal *T,*s,re,im;
20: PetscScalar *eigr,*eigi;
21: PetscInt i,n=10,l=2,k=5,ld;
22: PetscViewer viewer;
23: PetscBool verbose;
25: PetscFunctionBeginUser;
26: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
27: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
28: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GHIEP with compact storage - dimension %" PetscInt_FMT ".\n",n));
29: PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL));
30: PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
31: PetscCheck(l<=n && k<=n && l<=k,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Wrong value of dimensions");
32: PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
34: /* Create DS object */
35: PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
36: PetscCall(DSSetType(ds,DSGHIEP));
37: PetscCall(DSSetFromOptions(ds));
38: ld = n+2; /* test leading dimension larger than n */
39: PetscCall(DSAllocate(ds,ld));
40: PetscCall(DSSetDimensions(ds,n,l,k));
41: PetscCall(DSSetCompact(ds,PETSC_TRUE));
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: PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
50: /* Fill arrow-tridiagonal matrix */
51: PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T));
52: PetscCall(DSGetArrayReal(ds,DS_MAT_D,&s));
53: for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
54: for (i=k;i<n-1;i++) T[i+ld] = 1.0;
55: for (i=l;i<k;i++) T[i+2*ld] = 1.0;
56: T[2*ld+l+1] = -7; T[ld+k+1] = -7;
57: /* Signature matrix */
58: for (i=0;i<n;i++) s[i] = 1.0;
59: s[l+1] = -1.0;
60: s[k+1] = -1.0;
61: PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T));
62: PetscCall(DSRestoreArrayReal(ds,DS_MAT_D,&s));
63: if (l==0 && k==0) PetscCall(DSSetState(ds,DS_STATE_INTERMEDIATE));
64: else PetscCall(DSSetState(ds,DS_STATE_RAW));
65: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
66: PetscCall(DSView(ds,viewer));
68: /* Solve */
69: PetscCall(PetscCalloc2(n,&eigr,n,&eigi));
70: PetscCall(DSGetSlepcSC(ds,&sc));
71: sc->comparison = SlepcCompareLargestMagnitude;
72: sc->comparisonctx = NULL;
73: sc->map = NULL;
74: sc->mapobj = NULL;
75: PetscCall(DSSolve(ds,eigr,eigi));
76: PetscCall(DSSort(ds,eigr,eigi,NULL,NULL,NULL));
77: if (verbose) {
78: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
79: PetscCall(DSView(ds,viewer));
80: }
82: /* Print eigenvalues */
83: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n"));
84: for (i=0;i<n;i++) {
85: #if defined(PETSC_USE_COMPLEX)
86: re = PetscRealPart(eigr[i]);
87: im = PetscImaginaryPart(eigr[i]);
88: #else
89: re = eigr[i];
90: im = eigi[i];
91: #endif
92: if (PetscAbs(im)<1e-10) PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re));
93: else PetscCall(PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im));
94: }
95: PetscCall(PetscFree2(eigr,eigi));
96: PetscCall(DSDestroy(&ds));
97: PetscCall(SlepcFinalize());
98: return 0;
99: }
101: /*TEST
103: test:
104: suffix: 1
105: requires: !single
106: args: -ds_method {{0 1 2}}
107: filter: grep -v "solving the problem"
109: test:
110: suffix: 2
111: args: -n 5 -k 4 -l 4 -ds_method {{0 1 2}}
112: filter: grep -v "solving the problem"
114: TEST*/