Actual source code: test9.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 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;
25: DSGetLeadingDimension(ds,&ld);
26: DSGetDimensions(ds,&n,NULL,NULL,NULL);
27: DSGetArray(ds,mat,&X);
28: for (i=0;i<n;i++) {
29: aux = PetscAbsScalar(X[i+j*ld]);
30: nrm += aux*aux;
31: }
32: DSRestoreArray(ds,mat,&X);
33: *onrm = PetscSqrtReal(nrm);
34: PetscFunctionReturn(0);
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: SlepcInitialize(&argc,&argv,(char*)0,help);
49: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
50: PetscPrintf(PETSC_COMM_WORLD,"Solve a System of type GHEP - dimension %" PetscInt_FMT ".\n",n);
51: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
53: /* Create DS object */
54: DSCreate(PETSC_COMM_WORLD,&ds);
55: DSSetType(ds,DSGHEP);
56: DSSetFromOptions(ds);
57: ld = n+2; /* test leading dimension larger than n */
58: DSAllocate(ds,ld);
59: DSSetDimensions(ds,n,0,0);
61: /* Set up viewer */
62: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
63: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
64: DSView(ds,viewer);
65: PetscViewerPopFormat(viewer);
66: if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
68: /* Fill with a symmetric Toeplitz matrix */
69: DSGetArray(ds,DS_MAT_A,&A);
70: DSGetArray(ds,DS_MAT_B,&B);
71: for (i=0;i<n;i++) A[i+i*ld]=2.0;
72: for (j=1;j<3;j++) {
73: for (i=0;i<n-j;i++) { A[i+(i+j)*ld]=1.0; A[(i+j)+i*ld]=1.0; }
74: }
75: for (j=1;j<3;j++) { A[0+j*ld]=-1.0*(j+2); A[j+0*ld]=-1.0*(j+2); }
76: /* Diagonal matrix */
77: for (i=0;i<n;i++) B[i+i*ld]=0.1*(i+1);
78: DSRestoreArray(ds,DS_MAT_A,&A);
79: DSRestoreArray(ds,DS_MAT_B,&B);
80: DSSetState(ds,DS_STATE_RAW);
81: if (verbose) {
82: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
83: DSView(ds,viewer);
84: }
86: /* Solve */
87: PetscMalloc1(n,&eig);
88: PetscNew(&sc);
89: sc->comparison = SlepcCompareLargestMagnitude;
90: sc->comparisonctx = NULL;
91: sc->map = NULL;
92: sc->mapobj = NULL;
93: DSSetSlepcSC(ds,sc);
94: DSSolve(ds,eig,NULL);
95: DSSort(ds,eig,NULL,NULL,NULL,NULL);
96: if (verbose) {
97: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
98: DSView(ds,viewer);
99: }
101: /* Print eigenvalues */
102: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
103: for (i=0;i<n;i++) {
104: re = PetscRealPart(eig[i]);
105: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re);
106: }
108: /* Eigenvectors */
109: j = 0;
110: DSVectors(ds,DS_MAT_X,&j,NULL); /* all eigenvectors */
111: ComputeNorm(ds,DS_MAT_X,0,&nrm);
112: PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st vector = %.3f\n",(double)nrm);
113: DSVectors(ds,DS_MAT_X,NULL,NULL); /* all eigenvectors */
114: if (verbose) {
115: PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n");
116: DSView(ds,viewer);
117: }
119: PetscFree(eig);
120: PetscFree(sc);
121: DSDestroy(&ds);
122: SlepcFinalize();
123: return 0;
124: }
126: /*TEST
128: test:
129: suffix: 1
130: requires: !single
132: TEST*/