Actual source code: test24.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[] = "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: SlepcInitialize(&argc,&argv,(char*)0,help);
29: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
30: PetscOptionsGetInt(NULL,NULL,"-nini",&nini,NULL);
31: PetscOptionsGetInt(NULL,NULL,"-ncon",&ncon,NULL);
32: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem, n=%" PetscInt_FMT " nini=%" PetscInt_FMT " ncon=%" PetscInt_FMT "\n\n",n,nini,ncon);
34: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35: Compute the operator matrix that defines the eigensystem, Ax=kx
36: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
38: MatCreate(PETSC_COMM_WORLD,&A);
39: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
40: MatSetFromOptions(A);
41: MatSetUp(A);
43: MatGetOwnershipRange(A,&Istart,&Iend);
44: for (i=Istart;i<Iend;i++) {
45: if (i>0) MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
46: if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
47: MatSetValue(A,i,i,2.0,INSERT_VALUES);
48: }
49: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
50: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
52: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53: Create the eigensolver and set various options
54: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: EPSCreate(PETSC_COMM_WORLD,&eps);
56: EPSSetOperators(eps,A,NULL);
57: EPSSetProblemType(eps,EPS_HEP);
58: EPSSetType(eps,EPSLOBPCG);
59: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
60: EPSSetConvergenceTest(eps,EPS_CONV_ABS);
61: EPSSetDimensions(eps,nev,PETSC_DEFAULT,PETSC_DEFAULT);
62: EPSLOBPCGSetBlockSize(eps,nev);
63: EPSLOBPCGSetRestart(eps,0.7);
64: EPSSetTolerances(eps,1e-8,1200);
65: EPSSetFromOptions(eps);
67: MatCreateVecs(A,&t,NULL);
68: if (nini) {
69: VecDuplicateVecs(t,nini,&vi);
70: for (i=0;i<nini;i++) VecSetRandom(vi[i],NULL);
71: EPSSetInitialSpace(eps,nini,vi);
72: }
73: if (ncon) { /* constraints are exact eigenvectors of lowest eigenvalues */
74: alpha = PETSC_PI/(n+1);
75: beta = PetscSqrtReal(2.0/(n+1));
76: VecDuplicateVecs(t,ncon,&vc);
77: for (i=0;i<ncon;i++) {
78: for (j=0;j<n;j++) VecSetValue(vc[i],j,PetscSinReal(alpha*(j+1)*(i+1))*beta,INSERT_VALUES);
79: VecAssemblyBegin(vc[i]);
80: VecAssemblyEnd(vc[i]);
81: }
82: EPSSetDeflationSpace(eps,ncon,vc);
83: }
85: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86: Solve the eigensystem
87: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
89: EPSSolve(eps);
90: EPSGetType(eps,&type);
91: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n",type);
92: PetscObjectTypeCompare((PetscObject)eps,EPSLOBPCG,&flg);
93: if (flg) {
94: EPSLOBPCGGetLocking(eps,&lock);
95: if (lock) PetscPrintf(PETSC_COMM_WORLD," Using soft locking\n");
96: EPSLOBPCGGetRestart(eps,&restart);
97: PetscPrintf(PETSC_COMM_WORLD," LOBPCG Restart parameter=%.4g\n",(double)restart);
98: EPSLOBPCGGetBlockSize(eps,&bs);
99: PetscPrintf(PETSC_COMM_WORLD," LOBPCG Block size=%" PetscInt_FMT "\n",bs);
100: }
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Display solution and clean up
104: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
107: EPSDestroy(&eps);
108: MatDestroy(&A);
109: VecDestroyVecs(nini,&vi);
110: VecDestroyVecs(ncon,&vc);
111: VecDestroy(&t);
112: SlepcFinalize();
113: return 0;
114: }
116: /*TEST
118: testset:
119: args: -ncon 2
120: output_file: output/test24_1.out
121: test:
122: suffix: 1
123: requires: !single
124: test:
125: suffix: 2_cuda
126: args: -mat_type aijcusparse
127: requires: cuda !single
129: TEST*/