Actual source code: test9.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 DSGHEP.\n\n";
13: #include <slepcds.h>
15: /*
16: Compute the norm of the j-th column of matrix mat in ds
17: */
18: PetscErrorCode ComputeNorm(DS ds,DSMatType mat,PetscInt j,PetscReal *onrm)
19: {
20: PetscScalar *X;
21: PetscReal aux,nrm=0.0;
22: PetscInt i,n,ld;
24: PetscFunctionBeginUser;
25: PetscCall(DSGetLeadingDimension(ds,&ld));
26: PetscCall(DSGetDimensions(ds,&n,NULL,NULL,NULL));
27: PetscCall(DSGetArray(ds,mat,&X));
28: for (i=0;i<n;i++) {
29: aux = PetscAbsScalar(X[i+j*ld]);
30: nrm += aux*aux;
31: }
32: PetscCall(DSRestoreArray(ds,mat,&X));
33: *onrm = PetscSqrtReal(nrm);
34: PetscFunctionReturn(PETSC_SUCCESS);
35: }
37: int main(int argc,char **argv)
38: {
39: DS ds;
40: SlepcSC sc;
41: PetscReal re;
42: PetscScalar *A,*B,*eig;
43: PetscReal nrm;
44: PetscInt i,j,n=10,ld;
45: PetscViewer viewer;
46: PetscBool verbose;
48: PetscFunctionBeginUser;
49: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
50: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
51: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a System of type GHEP - dimension %" PetscInt_FMT ".\n",n));
52: PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
54: /* Create DS object */
55: PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
56: PetscCall(DSSetType(ds,DSGHEP));
57: PetscCall(DSSetFromOptions(ds));
58: ld = n+2; /* test leading dimension larger than n */
59: PetscCall(DSAllocate(ds,ld));
60: PetscCall(DSSetDimensions(ds,n,0,0));
62: /* Set up viewer */
63: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
64: PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
65: PetscCall(DSView(ds,viewer));
66: PetscCall(PetscViewerPopFormat(viewer));
67: if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
69: /* Fill with a symmetric Toeplitz matrix */
70: PetscCall(DSGetArray(ds,DS_MAT_A,&A));
71: PetscCall(DSGetArray(ds,DS_MAT_B,&B));
72: for (i=0;i<n;i++) A[i+i*ld]=2.0;
73: for (j=1;j<3;j++) {
74: for (i=0;i<n-j;i++) { A[i+(i+j)*ld]=1.0; A[(i+j)+i*ld]=1.0; }
75: }
76: for (j=1;j<3;j++) { A[0+j*ld]=-1.0*(j+2); A[j+0*ld]=-1.0*(j+2); }
77: /* Diagonal matrix */
78: for (i=0;i<n;i++) B[i+i*ld]=0.1*(i+1);
79: PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
80: PetscCall(DSRestoreArray(ds,DS_MAT_B,&B));
81: PetscCall(DSSetState(ds,DS_STATE_RAW));
82: if (verbose) {
83: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
84: PetscCall(DSView(ds,viewer));
85: }
87: /* Solve */
88: PetscCall(PetscMalloc1(n,&eig));
89: PetscCall(PetscNew(&sc));
90: sc->comparison = SlepcCompareLargestMagnitude;
91: sc->comparisonctx = NULL;
92: sc->map = NULL;
93: sc->mapobj = NULL;
94: PetscCall(DSSetSlepcSC(ds,sc));
95: PetscCall(DSSolve(ds,eig,NULL));
96: PetscCall(DSSort(ds,eig,NULL,NULL,NULL,NULL));
97: if (verbose) {
98: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
99: PetscCall(DSView(ds,viewer));
100: }
102: /* Print eigenvalues */
103: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n"));
104: for (i=0;i<n;i++) {
105: re = PetscRealPart(eig[i]);
106: PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re));
107: }
109: /* Eigenvectors */
110: j = 0;
111: PetscCall(DSVectors(ds,DS_MAT_X,&j,NULL)); /* all eigenvectors */
112: PetscCall(ComputeNorm(ds,DS_MAT_X,0,&nrm));
113: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st vector = %.3f\n",(double)nrm));
114: PetscCall(DSVectors(ds,DS_MAT_X,NULL,NULL)); /* all eigenvectors */
115: if (verbose) {
116: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n"));
117: PetscCall(DSView(ds,viewer));
118: }
120: PetscCall(PetscFree(eig));
121: PetscCall(PetscFree(sc));
122: PetscCall(DSDestroy(&ds));
123: PetscCall(SlepcFinalize());
124: return 0;
125: }
127: /*TEST
129: test:
130: suffix: 1
131: requires: !single
133: TEST*/