Actual source code: test13.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 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: SlepcInitialize(&argc,&argv,(char*)0,help);
25: PetscOptionsGetInt(NULL,NULL,"-maxbw",&maxbw,NULL);
26: PetscOptionsGetInt(NULL,NULL,"-nblks",&nblks,NULL);
27: n = maxbw*nblks;
28: bs = maxbw;
29: PetscPrintf(PETSC_COMM_WORLD,"Solve a block HEP Dense System - dimension %" PetscInt_FMT " (bandwidth=%" PetscInt_FMT ", blocks=%" PetscInt_FMT ").\n",n,maxbw,nblks);
30: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
32: /* Create DS object */
33: DSCreate(PETSC_COMM_WORLD,&ds);
34: DSSetType(ds,DSHEP);
35: DSSetMethod(ds,3); /* Select block divide-and-conquer */
36: DSSetBlockSize(ds,bs);
37: DSSetFromOptions(ds);
38: ld = n;
39: DSAllocate(ds,ld);
40: DSSetDimensions(ds,n,0,0);
42: /* Set up viewer */
43: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
44: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
45: DSView(ds,viewer);
46: PetscViewerPopFormat(viewer);
47: if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
49: /* Fill with a symmetric band Toeplitz matrix */
50: DSGetArray(ds,DS_MAT_A,&A);
51: for (i=0;i<n;i++) A[i+i*ld]=2.0;
52: for (j=1;j<=bs;j++) {
53: for (i=0;i<n-j;i++) { A[i+(i+j)*ld]=1.0; A[(i+j)+i*ld]=1.0; }
54: }
55: DSRestoreArray(ds,DS_MAT_A,&A);
56: DSSetState(ds,DS_STATE_RAW);
57: if (verbose) {
58: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
59: DSView(ds,viewer);
60: }
62: /* Solve */
63: PetscMalloc1(n,&eig);
64: DSGetSlepcSC(ds,&sc);
65: sc->comparison = SlepcCompareSmallestReal;
66: sc->comparisonctx = NULL;
67: sc->map = NULL;
68: sc->mapobj = NULL;
69: DSSolve(ds,eig,NULL);
70: DSSort(ds,eig,NULL,NULL,NULL,NULL);
71: if (verbose) {
72: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
73: DSView(ds,viewer);
74: }
76: /* Print eigenvalues */
77: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
78: for (i=0;i<n;i++) PetscViewerASCIIPrintf(viewer," %.5f\n",(double)PetscRealPart(eig[i]));
80: PetscFree(eig);
81: DSDestroy(&ds);
82: SlepcFinalize();
83: return 0;
84: }
86: /*TEST
88: test:
89: suffix: 1
90: requires: !complex !single
92: test:
93: suffix: 2
94: args: -nblks 1
95: requires: !complex
97: TEST*/