Actual source code: slepccontour.c
slepc-3.17.1 2022-04-11
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: #include <slepc/private/slepccontour.h>
12: #include <slepcblaslapack.h>
14: /*
15: SlepcContourDataCreate - Create a contour data structure.
17: Input Parameters:
18: n - the number of integration points
19: npart - number of partitions for the subcommunicator
20: parent - parent object
21: */
22: PetscErrorCode SlepcContourDataCreate(PetscInt n,PetscInt npart,PetscObject parent,SlepcContourData *contour)
23: {
24: PetscNew(contour);
25: (*contour)->parent = parent;
26: PetscSubcommCreate(PetscObjectComm(parent),&(*contour)->subcomm);
27: PetscSubcommSetNumber((*contour)->subcomm,npart);
28: PetscSubcommSetType((*contour)->subcomm,PETSC_SUBCOMM_INTERLACED);
29: PetscLogObjectMemory(parent,sizeof(PetscSubcomm));
30: (*contour)->npoints = n / npart;
31: if (n%npart > (*contour)->subcomm->color) (*contour)->npoints++;
32: PetscFunctionReturn(0);
33: }
35: /*
36: SlepcContourDataReset - Resets the KSP objects in a contour data structure,
37: and destroys any objects whose size depends on the problem size.
38: */
39: PetscErrorCode SlepcContourDataReset(SlepcContourData contour)
40: {
41: PetscInt i;
43: if (contour->ksp) {
44: for (i=0;i<contour->npoints;i++) KSPReset(contour->ksp[i]);
45: }
46: if (contour->pA) {
47: MatDestroyMatrices(contour->nmat,&contour->pA);
48: MatDestroyMatrices(contour->nmat,&contour->pP);
49: contour->nmat = 0;
50: }
51: VecScatterDestroy(&contour->scatterin);
52: VecDestroy(&contour->xsub);
53: VecDestroy(&contour->xdup);
54: PetscFunctionReturn(0);
55: }
57: /*
58: SlepcContourDataDestroy - Destroys the contour data structure.
59: */
60: PetscErrorCode SlepcContourDataDestroy(SlepcContourData *contour)
61: {
62: PetscInt i;
64: if (!(*contour)) PetscFunctionReturn(0);
65: if ((*contour)->ksp) {
66: for (i=0;i<(*contour)->npoints;i++) KSPDestroy(&(*contour)->ksp[i]);
67: PetscFree((*contour)->ksp);
68: }
69: PetscSubcommDestroy(&(*contour)->subcomm);
70: PetscFree((*contour));
71: *contour = NULL;
72: PetscFunctionReturn(0);
73: }
75: /*
76: SlepcContourRedundantMat - Creates redundant copies of the passed matrices in the subcomm.
78: Input Parameters:
79: nmat - the number of matrices
80: A - array of matrices
81: P - array of matrices (preconditioner)
82: */
83: PetscErrorCode SlepcContourRedundantMat(SlepcContourData contour,PetscInt nmat,Mat *A,Mat *P)
84: {
85: PetscInt i;
86: MPI_Comm child;
88: if (contour->pA) {
89: MatDestroyMatrices(contour->nmat,&contour->pA);
90: MatDestroyMatrices(contour->nmat,&contour->pP);
91: contour->nmat = 0;
92: }
93: if (contour->subcomm && contour->subcomm->n != 1) {
94: PetscSubcommGetChild(contour->subcomm,&child);
95: PetscCalloc1(nmat,&contour->pA);
96: for (i=0;i<nmat;i++) {
97: MatCreateRedundantMatrix(A[i],contour->subcomm->n,child,MAT_INITIAL_MATRIX,&contour->pA[i]);
98: PetscLogObjectParent(contour->parent,(PetscObject)contour->pA[i]);
99: }
100: if (P) {
101: PetscCalloc1(nmat,&contour->pP);
102: for (i=0;i<nmat;i++) {
103: MatCreateRedundantMatrix(P[i],contour->subcomm->n,child,MAT_INITIAL_MATRIX,&contour->pP[i]);
104: PetscLogObjectParent(contour->parent,(PetscObject)contour->pP[i]);
105: }
106: }
107: contour->nmat = nmat;
108: }
109: PetscFunctionReturn(0);
110: }
112: /*
113: SlepcContourScatterCreate - Creates a scatter context to communicate between a
114: regular vector and a vector xdup that can hold one duplicate per each subcommunicator
115: on the contiguous parent communicator. Also creates auxiliary vectors xdup and xsub
116: (the latter with the same layout as the redundant matrices in the subcommunicator).
118: Input Parameters:
119: v - the regular vector from which dimensions are taken
120: */
121: PetscErrorCode SlepcContourScatterCreate(SlepcContourData contour,Vec v)
122: {
123: IS is1,is2;
124: PetscInt i,j,k,m,mstart,mend,mlocal;
125: PetscInt *idx1,*idx2,mloc_sub;
126: MPI_Comm contpar,parent;
128: VecDestroy(&contour->xsub);
129: MatCreateVecsEmpty(contour->pA[0],&contour->xsub,NULL);
131: VecDestroy(&contour->xdup);
132: MatGetLocalSize(contour->pA[0],&mloc_sub,NULL);
133: PetscSubcommGetContiguousParent(contour->subcomm,&contpar);
134: VecCreate(contpar,&contour->xdup);
135: VecSetSizes(contour->xdup,mloc_sub,PETSC_DECIDE);
136: VecSetType(contour->xdup,((PetscObject)v)->type_name);
138: VecScatterDestroy(&contour->scatterin);
139: VecGetSize(v,&m);
140: VecGetOwnershipRange(v,&mstart,&mend);
141: mlocal = mend - mstart;
142: PetscMalloc2(contour->subcomm->n*mlocal,&idx1,contour->subcomm->n*mlocal,&idx2);
143: j = 0;
144: for (k=0;k<contour->subcomm->n;k++) {
145: for (i=mstart;i<mend;i++) {
146: idx1[j] = i;
147: idx2[j++] = i + m*k;
148: }
149: }
150: PetscSubcommGetParent(contour->subcomm,&parent);
151: ISCreateGeneral(parent,contour->subcomm->n*mlocal,idx1,PETSC_COPY_VALUES,&is1);
152: ISCreateGeneral(parent,contour->subcomm->n*mlocal,idx2,PETSC_COPY_VALUES,&is2);
153: VecScatterCreate(v,is1,contour->xdup,is2,&contour->scatterin);
154: ISDestroy(&is1);
155: ISDestroy(&is2);
156: PetscFree2(idx1,idx2);
157: PetscFunctionReturn(0);
158: }
160: /*
161: SlepcCISS_isGhost - Determine if any of the computed eigenpairs are spurious.
163: Input Parameters:
164: X - the matrix of eigenvectors (MATSEQDENSE)
165: n - the number of columns to consider
166: sigma - the singular values
167: thresh - threshold to decide whether a value is spurious
169: Output Parameter:
170: fl - array of n booleans
171: */
172: PetscErrorCode SlepcCISS_isGhost(Mat X,PetscInt n,PetscReal *sigma,PetscReal thresh,PetscBool *fl)
173: {
174: const PetscScalar *pX;
175: PetscInt i,j,m,ld;
176: PetscReal *tau,s1,s2,tau_max=0.0;
178: MatGetSize(X,&m,NULL);
179: MatDenseGetLDA(X,&ld);
180: PetscMalloc1(n,&tau);
181: MatDenseGetArrayRead(X,&pX);
182: for (j=0;j<n;j++) {
183: s1 = 0.0;
184: s2 = 0.0;
185: for (i=0;i<m;i++) {
186: s1 += PetscAbsScalar(PetscPowScalarInt(pX[i+j*ld],2));
187: s2 += PetscPowRealInt(PetscAbsScalar(pX[i+j*ld]),2)/sigma[i];
188: }
189: tau[j] = s1/s2;
190: tau_max = PetscMax(tau_max,tau[j]);
191: }
192: MatDenseRestoreArrayRead(X,&pX);
193: for (j=0;j<n;j++) fl[j] = (tau[j]>=thresh*tau_max)? PETSC_TRUE: PETSC_FALSE;
194: PetscFree(tau);
195: PetscFunctionReturn(0);
196: }
198: /*
199: SlepcCISS_BH_SVD - Compute SVD of block Hankel matrix and its rank.
201: Input Parameters:
202: H - block Hankel matrix obtained via CISS_BlockHankel()
203: ml - dimension of rows and columns, equal to M*L
204: delta - the tolerance used to determine the rank
206: Output Parameters:
207: sigma - computed singular values
208: rank - the rank of H
209: */
210: PetscErrorCode SlepcCISS_BH_SVD(PetscScalar *H,PetscInt ml,PetscReal delta,PetscReal *sigma,PetscInt *rank)
211: {
212: PetscInt i;
213: PetscBLASInt m,n,lda,ldu,ldvt,lwork,info;
214: PetscScalar *work;
215: #if defined(PETSC_USE_COMPLEX)
216: PetscReal *rwork;
217: #endif
219: PetscMalloc1(5*ml,&work);
220: #if defined(PETSC_USE_COMPLEX)
221: PetscMalloc1(5*ml,&rwork);
222: #endif
223: PetscBLASIntCast(ml,&m);
224: n = m; lda = m; ldu = m; ldvt = m; lwork = 5*m;
225: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
226: #if defined(PETSC_USE_COMPLEX)
227: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,H,&lda,sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,rwork,&info));
228: #else
229: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,H,&lda,sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,&info));
230: #endif
231: SlepcCheckLapackInfo("gesvd",info);
232: PetscFPTrapPop();
233: (*rank) = 0;
234: for (i=0;i<ml;i++) {
235: if (sigma[i]/PetscMax(sigma[0],1.0)>delta) (*rank)++;
236: }
237: PetscFree(work);
238: #if defined(PETSC_USE_COMPLEX)
239: PetscFree(rwork);
240: #endif
241: PetscFunctionReturn(0);
242: }