Actual source code: slepccontour.c
slepc-3.16.0 2021-09-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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: {
27: PetscNew(contour);
28: (*contour)->parent = parent;
29: PetscSubcommCreate(PetscObjectComm(parent),&(*contour)->subcomm);
30: PetscSubcommSetNumber((*contour)->subcomm,npart);
31: PetscSubcommSetType((*contour)->subcomm,PETSC_SUBCOMM_INTERLACED);
32: PetscLogObjectMemory(parent,sizeof(PetscSubcomm));
33: (*contour)->npoints = n / npart;
34: if (n%npart > (*contour)->subcomm->color) (*contour)->npoints++;
35: return(0);
36: }
38: /*
39: SlepcContourDataReset - Resets the KSP objects in a contour data structure,
40: and destroys any objects whose size depends on the problem size.
41: */
42: PetscErrorCode SlepcContourDataReset(SlepcContourData contour)
43: {
45: PetscInt i;
48: if (contour->ksp) {
49: for (i=0;i<contour->npoints;i++) {
50: KSPReset(contour->ksp[i]);
51: }
52: }
53: if (contour->pA) {
54: for (i=0;i<contour->nmat;i++) {
55: MatDestroy(&contour->pA[i]);
56: }
57: PetscFree(contour->pA);
58: contour->pA = NULL;
59: contour->nmat = 0;
60: }
61: VecScatterDestroy(&contour->scatterin);
62: VecDestroy(&contour->xsub);
63: VecDestroy(&contour->xdup);
64: return(0);
65: }
67: /*
68: SlepcContourDataDestroy - Destroys the contour data structure.
69: */
70: PetscErrorCode SlepcContourDataDestroy(SlepcContourData *contour)
71: {
73: PetscInt i;
76: if (!(*contour)) return(0);
77: if ((*contour)->ksp) {
78: for (i=0;i<(*contour)->npoints;i++) {
79: KSPDestroy(&(*contour)->ksp[i]);
80: }
81: PetscFree((*contour)->ksp);
82: }
83: PetscSubcommDestroy(&(*contour)->subcomm);
84: PetscFree((*contour));
85: *contour = NULL;
86: return(0);
87: }
89: /*
90: SlepcContourRedundantMat - Creates redundant copies of the passed matrices in the subcomm.
92: Input Parameters:
93: nmat - the number of matrices
94: A - array of matrices
95: */
96: PetscErrorCode SlepcContourRedundantMat(SlepcContourData contour,PetscInt nmat,Mat *A)
97: {
99: PetscInt i;
102: if (contour->pA) {
103: for (i=0;i<contour->nmat;i++) {
104: MatDestroy(&contour->pA[i]);
105: }
106: PetscFree(contour->pA);
107: contour->pA = NULL;
108: contour->nmat = 0;
109: }
110: if (contour->subcomm && contour->subcomm->n != 1) {
111: PetscCalloc1(nmat,&contour->pA);
112: for (i=0;i<nmat;i++) {
113: MatCreateRedundantMatrix(A[i],contour->subcomm->n,PetscSubcommChild(contour->subcomm),MAT_INITIAL_MATRIX,&contour->pA[i]);
114: PetscLogObjectParent(contour->parent,(PetscObject)contour->pA[i]);
115: }
116: contour->nmat = nmat;
117: }
118: return(0);
119: }
121: /*
122: SlepcContourScatterCreate - Creates a scatter context to communicate between a
123: regular vector and a vector xdup that can hold one duplicate per each subcommunicator
124: on the contiguous parent communicator. Also creates auxiliary vectors xdup and xsub
125: (the latter with the same layout as the redundant matrices in the subcommunicator).
127: Input Parameters:
128: v - the regular vector from which dimensions are taken
129: */
130: PetscErrorCode SlepcContourScatterCreate(SlepcContourData contour,Vec v)
131: {
133: IS is1,is2;
134: PetscInt i,j,k,m,mstart,mend,mlocal;
135: PetscInt *idx1,*idx2,mloc_sub;
138: VecDestroy(&contour->xsub);
139: MatCreateVecsEmpty(contour->pA[0],&contour->xsub,NULL);
141: VecDestroy(&contour->xdup);
142: MatGetLocalSize(contour->pA[0],&mloc_sub,NULL);
143: VecCreate(PetscSubcommContiguousParent(contour->subcomm),&contour->xdup);
144: VecSetSizes(contour->xdup,mloc_sub,PETSC_DECIDE);
145: VecSetType(contour->xdup,((PetscObject)v)->type_name);
147: VecScatterDestroy(&contour->scatterin);
148: VecGetSize(v,&m);
149: VecGetOwnershipRange(v,&mstart,&mend);
150: mlocal = mend - mstart;
151: PetscMalloc2(contour->subcomm->n*mlocal,&idx1,contour->subcomm->n*mlocal,&idx2);
152: j = 0;
153: for (k=0;k<contour->subcomm->n;k++) {
154: for (i=mstart;i<mend;i++) {
155: idx1[j] = i;
156: idx2[j++] = i + m*k;
157: }
158: }
159: ISCreateGeneral(PetscSubcommParent(contour->subcomm),contour->subcomm->n*mlocal,idx1,PETSC_COPY_VALUES,&is1);
160: ISCreateGeneral(PetscSubcommParent(contour->subcomm),contour->subcomm->n*mlocal,idx2,PETSC_COPY_VALUES,&is2);
161: VecScatterCreate(v,is1,contour->xdup,is2,&contour->scatterin);
162: ISDestroy(&is1);
163: ISDestroy(&is2);
164: PetscFree2(idx1,idx2);
165: return(0);
166: }
168: /*
169: SlepcCISS_isGhost - Determine if any of the computed eigenpairs are spurious.
171: Input Parameters:
172: X - the matrix of eigenvectors (MATSEQDENSE)
173: n - the number of columns to consider
174: sigma - the singular values
175: thresh - threshold to decide whether a value is spurious
177: Output Parameter:
178: fl - array of n booleans
179: */
180: PetscErrorCode SlepcCISS_isGhost(Mat X,PetscInt n,PetscReal *sigma,PetscReal thresh,PetscBool *fl)
181: {
182: PetscErrorCode ierr;
183: const PetscScalar *pX;
184: PetscInt i,j,m,ld;
185: PetscReal *tau,s1,s2,tau_max=0.0;
188: MatGetSize(X,&m,NULL);
189: MatDenseGetLDA(X,&ld);
190: PetscMalloc1(n,&tau);
191: MatDenseGetArrayRead(X,&pX);
192: for (j=0;j<n;j++) {
193: s1 = 0.0;
194: s2 = 0.0;
195: for (i=0;i<m;i++) {
196: s1 += PetscAbsScalar(PetscPowScalarInt(pX[i+j*ld],2));
197: s2 += PetscPowRealInt(PetscAbsScalar(pX[i+j*ld]),2)/sigma[i];
198: }
199: tau[j] = s1/s2;
200: tau_max = PetscMax(tau_max,tau[j]);
201: }
202: MatDenseRestoreArrayRead(X,&pX);
203: for (j=0;j<n;j++) fl[j] = (tau[j]>=thresh*tau_max)? PETSC_TRUE: PETSC_FALSE;
204: PetscFree(tau);
205: return(0);
206: }
208: /*
209: SlepcCISS_BH_SVD - Compute SVD of block Hankel matrix and its rank.
211: Input Parameters:
212: H - block Hankel matrix obtained via CISS_BlockHankel()
213: ml - dimension of rows and columns, equal to M*L
214: delta - the tolerance used to determine the rank
216: Output Parameters:
217: sigma - computed singular values
218: rank - the rank of H
219: */
220: PetscErrorCode SlepcCISS_BH_SVD(PetscScalar *H,PetscInt ml,PetscReal delta,PetscReal *sigma,PetscInt *rank)
221: {
223: PetscInt i;
224: PetscBLASInt m,n,lda,ldu,ldvt,lwork,info;
225: PetscScalar *work;
226: #if defined(PETSC_USE_COMPLEX)
227: PetscReal *rwork;
228: #endif
231: PetscMalloc1(5*ml,&work);
232: #if defined(PETSC_USE_COMPLEX)
233: PetscMalloc1(5*ml,&rwork);
234: #endif
235: PetscBLASIntCast(ml,&m);
236: n = m; lda = m; ldu = m; ldvt = m; lwork = 5*m;
237: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
238: #if defined(PETSC_USE_COMPLEX)
239: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,H,&lda,sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,rwork,&info));
240: #else
241: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&m,&n,H,&lda,sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,&info));
242: #endif
243: SlepcCheckLapackInfo("gesvd",info);
244: PetscFPTrapPop();
245: (*rank) = 0;
246: for (i=0;i<ml;i++) {
247: if (sigma[i]/PetscMax(sigma[0],1.0)>delta) (*rank)++;
248: }
249: PetscFree(work);
250: #if defined(PETSC_USE_COMPLEX)
251: PetscFree(rwork);
252: #endif
253: return(0);
254: }