Actual source code: test24.c
slepc-3.20.2 2024-03-15
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[] = "Eigenproblem for the 1-D Laplacian with constraints. "
12: "Based on ex1.\n"
13: "The command line options are:\n"
14: " -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n\n";
16: #include <slepceps.h>
18: int main(int argc,char **argv)
19: {
20: Mat A;
21: EPS eps;
22: EPSType type;
23: Vec *vi=NULL,*vc=NULL,t;
24: PetscInt n=30,nev=4,i,j,Istart,Iend,nini=0,ncon=0,bs;
25: PetscReal alpha,beta,restart;
26: PetscBool flg,lock;
28: PetscFunctionBeginUser;
29: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
30: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
31: PetscCall(PetscOptionsGetInt(NULL,NULL,"-nini",&nini,NULL));
32: PetscCall(PetscOptionsGetInt(NULL,NULL,"-ncon",&ncon,NULL));
33: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem, n=%" PetscInt_FMT " nini=%" PetscInt_FMT " ncon=%" PetscInt_FMT "\n\n",n,nini,ncon));
35: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: Compute the operator matrix that defines the eigensystem, Ax=kx
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
40: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
41: PetscCall(MatSetFromOptions(A));
42: PetscCall(MatSetUp(A));
44: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
45: for (i=Istart;i<Iend;i++) {
46: if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
47: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
48: PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
49: }
50: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
51: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Create the eigensolver and set various options
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
57: PetscCall(EPSSetOperators(eps,A,NULL));
58: PetscCall(EPSSetProblemType(eps,EPS_HEP));
59: PetscCall(EPSSetType(eps,EPSLOBPCG));
60: PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
61: PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_ABS));
62: PetscCall(EPSSetDimensions(eps,nev,PETSC_DEFAULT,PETSC_DEFAULT));
63: PetscCall(EPSLOBPCGSetBlockSize(eps,nev));
64: PetscCall(EPSLOBPCGSetRestart(eps,0.7));
65: PetscCall(EPSSetTolerances(eps,1e-8,1200));
66: PetscCall(EPSSetFromOptions(eps));
68: PetscCall(MatCreateVecs(A,&t,NULL));
69: if (nini) {
70: PetscCall(VecDuplicateVecs(t,nini,&vi));
71: for (i=0;i<nini;i++) PetscCall(VecSetRandom(vi[i],NULL));
72: PetscCall(EPSSetInitialSpace(eps,nini,vi));
73: }
74: if (ncon) { /* constraints are exact eigenvectors of lowest eigenvalues */
75: alpha = PETSC_PI/(n+1);
76: beta = PetscSqrtReal(2.0/(n+1));
77: PetscCall(VecDuplicateVecs(t,ncon,&vc));
78: for (i=0;i<ncon;i++) {
79: for (j=0;j<n;j++) PetscCall(VecSetValue(vc[i],j,PetscSinReal(alpha*(j+1)*(i+1))*beta,INSERT_VALUES));
80: PetscCall(VecAssemblyBegin(vc[i]));
81: PetscCall(VecAssemblyEnd(vc[i]));
82: }
83: PetscCall(EPSSetDeflationSpace(eps,ncon,vc));
84: }
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Solve the eigensystem
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: PetscCall(EPSSolve(eps));
91: PetscCall(EPSGetType(eps,&type));
92: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n",type));
93: PetscCall(PetscObjectTypeCompare((PetscObject)eps,EPSLOBPCG,&flg));
94: if (flg) {
95: PetscCall(EPSLOBPCGGetLocking(eps,&lock));
96: if (lock) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Using soft locking\n"));
97: PetscCall(EPSLOBPCGGetRestart(eps,&restart));
98: PetscCall(PetscPrintf(PETSC_COMM_WORLD," LOBPCG Restart parameter=%.4g\n",(double)restart));
99: PetscCall(EPSLOBPCGGetBlockSize(eps,&bs));
100: PetscCall(PetscPrintf(PETSC_COMM_WORLD," LOBPCG Block size=%" PetscInt_FMT "\n",bs));
101: }
103: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: Display solution and clean up
105: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
108: PetscCall(EPSDestroy(&eps));
109: PetscCall(MatDestroy(&A));
110: PetscCall(VecDestroyVecs(nini,&vi));
111: PetscCall(VecDestroyVecs(ncon,&vc));
112: PetscCall(VecDestroy(&t));
113: PetscCall(SlepcFinalize());
114: return 0;
115: }
117: /*TEST
119: testset:
120: args: -ncon 2
121: output_file: output/test24_1.out
122: test:
123: suffix: 1
124: requires: !single
125: test:
126: suffix: 2_cuda
127: args: -mat_type aijcusparse
128: requires: cuda !single
130: TEST*/