Actual source code: test24.c
slepc-3.18.2 2023-01-26
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;
29: SlepcInitialize(&argc,&argv,(char*)0,help);
30: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
31: PetscOptionsGetInt(NULL,NULL,"-nini",&nini,NULL);
32: PetscOptionsGetInt(NULL,NULL,"-ncon",&ncon,NULL);
33: 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: MatCreate(PETSC_COMM_WORLD,&A);
40: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
41: MatSetFromOptions(A);
42: MatSetUp(A);
44: MatGetOwnershipRange(A,&Istart,&Iend);
45: for (i=Istart;i<Iend;i++) {
46: if (i>0) MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
47: if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
48: MatSetValue(A,i,i,2.0,INSERT_VALUES);
49: }
50: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
51: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Create the eigensolver and set various options
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: EPSCreate(PETSC_COMM_WORLD,&eps);
57: EPSSetOperators(eps,A,NULL);
58: EPSSetProblemType(eps,EPS_HEP);
59: EPSSetType(eps,EPSLOBPCG);
60: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
61: EPSSetConvergenceTest(eps,EPS_CONV_ABS);
62: EPSSetDimensions(eps,nev,PETSC_DEFAULT,PETSC_DEFAULT);
63: EPSLOBPCGSetBlockSize(eps,nev);
64: EPSLOBPCGSetRestart(eps,0.7);
65: EPSSetTolerances(eps,1e-8,1200);
66: EPSSetFromOptions(eps);
68: MatCreateVecs(A,&t,NULL);
69: if (nini) {
70: VecDuplicateVecs(t,nini,&vi);
71: for (i=0;i<nini;i++) VecSetRandom(vi[i],NULL);
72: 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: VecDuplicateVecs(t,ncon,&vc);
78: for (i=0;i<ncon;i++) {
79: for (j=0;j<n;j++) VecSetValue(vc[i],j,PetscSinReal(alpha*(j+1)*(i+1))*beta,INSERT_VALUES);
80: VecAssemblyBegin(vc[i]);
81: VecAssemblyEnd(vc[i]);
82: }
83: EPSSetDeflationSpace(eps,ncon,vc);
84: }
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Solve the eigensystem
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: EPSSolve(eps);
91: EPSGetType(eps,&type);
92: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n",type);
93: PetscObjectTypeCompare((PetscObject)eps,EPSLOBPCG,&flg);
94: if (flg) {
95: EPSLOBPCGGetLocking(eps,&lock);
96: if (lock) PetscPrintf(PETSC_COMM_WORLD," Using soft locking\n");
97: EPSLOBPCGGetRestart(eps,&restart);
98: PetscPrintf(PETSC_COMM_WORLD," LOBPCG Restart parameter=%.4g\n",(double)restart);
99: EPSLOBPCGGetBlockSize(eps,&bs);
100: PetscPrintf(PETSC_COMM_WORLD," LOBPCG Block size=%" PetscInt_FMT "\n",bs);
101: }
103: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: Display solution and clean up
105: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
108: EPSDestroy(&eps);
109: MatDestroy(&A);
110: VecDestroyVecs(nini,&vi);
111: VecDestroyVecs(ncon,&vc);
112: VecDestroy(&t);
113: 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*/