Actual source code: test13.c
slepc-3.21.1 2024-04-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[] = "Solve a quadratic problem with CISS.\n\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
14: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
16: #include <slepcpep.h>
18: int main(int argc,char **argv)
19: {
20: Mat M,C,K,A[3];
21: PEP pep;
22: RG rg;
23: KSP *ksp;
24: PC pc;
25: PEPCISSExtraction ext;
26: PetscInt N,n=10,m,Istart,Iend,II,i,j,nsolve;
27: PetscBool flg;
29: PetscFunctionBeginUser;
30: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
31: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
32: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flg));
33: if (!flg) m=n;
34: N = n*m;
35: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
37: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
39: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
41: /* K is the 2-D Laplacian */
42: PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
43: PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N));
44: PetscCall(MatSetFromOptions(K));
45: PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
46: for (II=Istart;II<Iend;II++) {
47: i = II/n; j = II-i*n;
48: if (i>0) PetscCall(MatSetValue(K,II,II-n,-1.0,INSERT_VALUES));
49: if (i<m-1) PetscCall(MatSetValue(K,II,II+n,-1.0,INSERT_VALUES));
50: if (j>0) PetscCall(MatSetValue(K,II,II-1,-1.0,INSERT_VALUES));
51: if (j<n-1) PetscCall(MatSetValue(K,II,II+1,-1.0,INSERT_VALUES));
52: PetscCall(MatSetValue(K,II,II,4.0,INSERT_VALUES));
53: }
54: PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
55: PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
57: /* C is the 1-D Laplacian on horizontal lines */
58: PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
59: PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N));
60: PetscCall(MatSetFromOptions(C));
61: PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
62: for (II=Istart;II<Iend;II++) {
63: i = II/n; j = II-i*n;
64: if (j>0) PetscCall(MatSetValue(C,II,II-1,-1.0,INSERT_VALUES));
65: if (j<n-1) PetscCall(MatSetValue(C,II,II+1,-1.0,INSERT_VALUES));
66: PetscCall(MatSetValue(C,II,II,2.0,INSERT_VALUES));
67: }
68: PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
69: PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
71: /* M is a diagonal matrix */
72: PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
73: PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N));
74: PetscCall(MatSetFromOptions(M));
75: PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
76: for (II=Istart;II<Iend;II++) PetscCall(MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES));
77: PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
78: PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
80: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: Create the eigensolver and solve the eigensystem
82: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
85: A[0] = K; A[1] = C; A[2] = M;
86: PetscCall(PEPSetOperators(pep,3,A));
87: PetscCall(PEPSetProblemType(pep,PEP_GENERAL));
89: /* customize polynomial eigensolver; set runtime options */
90: PetscCall(PEPSetType(pep,PEPCISS));
91: PetscCall(PEPGetRG(pep,&rg));
92: PetscCall(RGSetType(rg,RGELLIPSE));
93: PetscCall(RGEllipseSetParameters(rg,PetscCMPLX(-0.1,0.3),0.1,0.25));
94: PetscCall(PEPCISSSetSizes(pep,24,PETSC_DEFAULT,PETSC_DEFAULT,1,PETSC_DEFAULT,PETSC_TRUE));
95: PetscCall(PEPCISSGetKSPs(pep,&nsolve,&ksp));
96: for (i=0;i<nsolve;i++) {
97: PetscCall(KSPSetTolerances(ksp[i],1e-12,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
98: PetscCall(KSPSetType(ksp[i],KSPPREONLY));
99: PetscCall(KSPGetPC(ksp[i],&pc));
100: PetscCall(PCSetType(pc,PCLU));
101: }
102: PetscCall(PEPSetFromOptions(pep));
104: /* solve */
105: PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPCISS,&flg));
106: if (flg) {
107: PetscCall(PEPCISSGetExtraction(pep,&ext));
108: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Running CISS with %" PetscInt_FMT " KSP solvers (%s extraction)\n",nsolve,PEPCISSExtractions[ext]));
109: }
110: PetscCall(PEPSolve(pep));
112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: Display solution and clean up
114: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116: PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
117: PetscCall(PEPDestroy(&pep));
118: PetscCall(MatDestroy(&M));
119: PetscCall(MatDestroy(&C));
120: PetscCall(MatDestroy(&K));
121: PetscCall(SlepcFinalize());
122: return 0;
123: }
125: /*TEST
127: build:
128: requires: complex
130: test:
131: suffix: 1
132: requires: complex
134: TEST*/