Actual source code: test13.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 DSHEP with block size larger than one.\n\n";
13: #include <slepcds.h>
15: int main(int argc,char **argv)
16: {
17: DS ds;
18: SlepcSC sc;
19: PetscScalar *A,*eig;
20: PetscInt i,j,n,ld,bs,maxbw=3,nblks=8;
21: PetscViewer viewer;
22: PetscBool verbose;
24: PetscFunctionBeginUser;
25: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
26: PetscCall(PetscOptionsGetInt(NULL,NULL,"-maxbw",&maxbw,NULL));
27: PetscCall(PetscOptionsGetInt(NULL,NULL,"-nblks",&nblks,NULL));
28: n = maxbw*nblks;
29: bs = maxbw;
30: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a block HEP Dense System - dimension %" PetscInt_FMT " (bandwidth=%" PetscInt_FMT ", blocks=%" PetscInt_FMT ").\n",n,maxbw,nblks));
31: PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
33: /* Create DS object */
34: PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
35: PetscCall(DSSetType(ds,DSHEP));
36: PetscCall(DSSetMethod(ds,3)); /* Select block divide-and-conquer */
37: PetscCall(DSSetBlockSize(ds,bs));
38: PetscCall(DSSetFromOptions(ds));
39: ld = n;
40: PetscCall(DSAllocate(ds,ld));
41: PetscCall(DSSetDimensions(ds,n,0,0));
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 a symmetric band Toeplitz matrix */
51: PetscCall(DSGetArray(ds,DS_MAT_A,&A));
52: for (i=0;i<n;i++) A[i+i*ld]=2.0;
53: for (j=1;j<=bs;j++) {
54: for (i=0;i<n-j;i++) { A[i+(i+j)*ld]=1.0; A[(i+j)+i*ld]=1.0; }
55: }
56: PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
57: PetscCall(DSSetState(ds,DS_STATE_RAW));
58: if (verbose) {
59: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
60: PetscCall(DSView(ds,viewer));
61: }
63: /* Solve */
64: PetscCall(PetscMalloc1(n,&eig));
65: PetscCall(DSGetSlepcSC(ds,&sc));
66: sc->comparison = SlepcCompareSmallestReal;
67: sc->comparisonctx = NULL;
68: sc->map = NULL;
69: sc->mapobj = NULL;
70: PetscCall(DSSolve(ds,eig,NULL));
71: PetscCall(DSSort(ds,eig,NULL,NULL,NULL,NULL));
72: if (verbose) {
73: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
74: PetscCall(DSView(ds,viewer));
75: }
77: /* Print eigenvalues */
78: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n"));
79: for (i=0;i<n;i++) PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)PetscRealPart(eig[i])));
81: PetscCall(PetscFree(eig));
82: PetscCall(DSDestroy(&ds));
83: PetscCall(SlepcFinalize());
84: return 0;
85: }
87: /*TEST
89: test:
90: suffix: 1
91: requires: !complex !single
93: test:
94: suffix: 2
95: args: -nblks 1
96: requires: !complex
98: TEST*/