Actual source code: dspriv.c
1: /*
2: Private DS routines
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc-private/dsimpl.h> /*I "slepcds.h" I*/
25: #include <slepcblaslapack.h>
27: PetscErrorCode SlepcDenseMatProd(PetscScalar *C, PetscInt _ldC, PetscScalar b, PetscScalar a, const PetscScalar *A, PetscInt _ldA, PetscInt rA, PetscInt cA, PetscBool At, const PetscScalar *B, PetscInt _ldB, PetscInt rB, PetscInt cB, PetscBool Bt);
31: PetscErrorCode DSAllocateMat_Private(DS ds,DSMatType m)
32: {
33: PetscInt sz;
37: if (m==DS_MAT_T) sz = 3*ds->ld*sizeof(PetscScalar);
38: else if (m==DS_MAT_D) sz = ds->ld*sizeof(PetscScalar);
39: else sz = ds->ld*ds->ld*sizeof(PetscScalar);
40: if (ds->mat[m]) { PetscFree(ds->mat[m]); }
41: else { PetscLogObjectMemory(ds,sz); }
42: PetscMalloc(sz,&ds->mat[m]);
43: PetscMemzero(ds->mat[m],sz);
44: return(0);
45: }
49: PetscErrorCode DSAllocateMatReal_Private(DS ds,DSMatType m)
50: {
51: PetscInt sz;
55: if (m==DS_MAT_T) sz = 3*ds->ld*sizeof(PetscReal);
56: else if (m==DS_MAT_D) sz = ds->ld*sizeof(PetscReal);
57: else sz = ds->ld*ds->ld*sizeof(PetscReal);
58: if (!ds->rmat[m]) {
59: PetscLogObjectMemory(ds,sz);
60: PetscMalloc(sz,&ds->rmat[m]);
61: }
62: PetscMemzero(ds->rmat[m],sz);
63: return(0);
64: }
68: PetscErrorCode DSAllocateWork_Private(DS ds,PetscInt s,PetscInt r,PetscInt i)
69: {
73: if (s>ds->lwork) {
74: PetscFree(ds->work);
75: PetscMalloc(s*sizeof(PetscScalar),&ds->work);
76: PetscLogObjectMemory(ds,(s-ds->lwork)*sizeof(PetscScalar));
77: ds->lwork = s;
78: }
79: if (r>ds->lrwork) {
80: PetscFree(ds->rwork);
81: PetscMalloc(r*sizeof(PetscReal),&ds->rwork);
82: PetscLogObjectMemory(ds,(r-ds->lrwork)*sizeof(PetscReal));
83: ds->lrwork = r;
84: }
85: if (i>ds->liwork) {
86: PetscFree(ds->iwork);
87: PetscMalloc(i*sizeof(PetscBLASInt),&ds->iwork);
88: PetscLogObjectMemory(ds,(i-ds->liwork)*sizeof(PetscBLASInt));
89: ds->liwork = i;
90: }
91: return(0);
92: }
96: PetscErrorCode DSViewMat_Private(DS ds,PetscViewer viewer,DSMatType m)
97: {
98: PetscErrorCode ierr;
99: PetscInt i,j,rows,cols;
100: PetscScalar *v;
101: PetscViewerFormat format;
102: #if defined(PETSC_USE_COMPLEX)
103: PetscBool allreal = PETSC_TRUE;
104: #endif
107: PetscViewerGetFormat(viewer,&format);
108: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
109: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
110: if (ds->state==DS_STATE_TRUNCATED && m>=DS_MAT_Q) rows = ds->t;
111: else rows = (m==DS_MAT_A && ds->extrarow)? ds->n+1: ds->n;
112: cols = (ds->m!=0)? ds->m: ds->n;
113: #if defined(PETSC_USE_COMPLEX)
114: /* determine if matrix has all real values */
115: v = ds->mat[m];
116: for (i=0;i<rows;i++)
117: for (j=0;j<cols;j++)
118: if (PetscImaginaryPart(v[i+j*ds->ld])) { allreal = PETSC_FALSE; break; }
119: #endif
120: if (format == PETSC_VIEWER_ASCII_MATLAB) {
121: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",rows,cols);
122: PetscViewerASCIIPrintf(viewer,"%s = [\n",DSMatName[m]);
123: } else {
124: PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[m]);
125: }
127: for (i=0;i<rows;i++) {
128: v = ds->mat[m]+i;
129: for (j=0;j<cols;j++) {
130: #if defined(PETSC_USE_COMPLEX)
131: if (allreal) {
132: PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));
133: } else {
134: PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei ",PetscRealPart(*v),PetscImaginaryPart(*v));
135: }
136: #else
137: PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);
138: #endif
139: v += ds->ld;
140: }
141: PetscViewerASCIIPrintf(viewer,"\n");
142: }
144: if (format == PETSC_VIEWER_ASCII_MATLAB) {
145: PetscViewerASCIIPrintf(viewer,"];\n");
146: }
147: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
148: PetscViewerFlush(viewer);
149: return(0);
150: }
154: PetscErrorCode DSSortEigenvalues_Private(DS ds,PetscScalar *wr,PetscScalar *wi,PetscInt *perm,PetscBool isghiep)
155: {
157: PetscScalar re,im,wi0;
158: PetscInt i,j,result,tmp1,tmp2=0,d=1;
161: for (i=0;i<ds->n;i++) perm[i] = i;
162: /* insertion sort */
163: i=ds->l+1;
164: #if !defined(PETSC_USE_COMPLEX)
165: if (wi && wi[perm[i-1]]!=0.0) i++; /* initial value is complex */
166: #else
167: if (isghiep && PetscImaginaryPart(wr[perm[i-1]])!=0.0) i++;
168: #endif
169: for (;i<ds->n;i+=d) {
170: re = wr[perm[i]];
171: if (wi) im = wi[perm[i]];
172: else im = 0.0;
173: tmp1 = perm[i];
174: #if !defined(PETSC_USE_COMPLEX)
175: if (im!=0.0) { d = 2; tmp2 = perm[i+1]; }
176: else d = 1;
177: #else
178: if (isghiep && PetscImaginaryPart(re)!=0.0) { d = 2; tmp2 = perm[i+1]; }
179: else d = 1;
180: #endif
181: j = i-1;
182: if (wi) wi0 = wi[perm[j]];
183: else wi0 = 0.0;
184: (*ds->comp_fun)(re,im,wr[perm[j]],wi0,&result,ds->comp_ctx);
185: while (result<0 && j>=ds->l) {
186: perm[j+d] = perm[j];
187: j--;
188: #if !defined(PETSC_USE_COMPLEX)
189: if (wi && wi[perm[j+1]]!=0)
190: #else
191: if (isghiep && PetscImaginaryPart(wr[perm[j+1]])!=0)
192: #endif
193: { perm[j+d] = perm[j]; j--; }
195: if (j>=ds->l) {
196: if (wi) wi0 = wi[perm[j]];
197: else wi0 = 0.0;
198: (*ds->comp_fun)(re,im,wr[perm[j]],wi0,&result,ds->comp_ctx);
199: }
200: }
201: perm[j+1] = tmp1;
202: if(d==2) perm[j+2] = tmp2;
203: }
204: return(0);
205: }
209: PetscErrorCode DSSortEigenvaluesReal_Private(DS ds,PetscReal *eig,PetscInt *perm)
210: {
212: PetscScalar re;
213: PetscInt i,j,result,tmp,l,n;
216: n = ds->n;
217: l = ds->l;
218: for (i=0;i<n;i++) perm[i] = i;
219: /* insertion sort */
220: for (i=l+1;i<n;i++) {
221: re = eig[perm[i]];
222: j = i-1;
223: (*ds->comp_fun)(re,0.0,eig[perm[j]],0.0,&result,ds->comp_ctx);
224: while (result<0 && j>=l) {
225: tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
226: if (j>=l) {
227: (*ds->comp_fun)(re,0.0,eig[perm[j]],0.0,&result,ds->comp_ctx);
228: }
229: }
230: }
231: return(0);
232: }
236: /*
237: DSCopyMatrix_Private - Copies the trailing block of a matrix (from
238: rows/columns l to n).
239: */
240: PetscErrorCode DSCopyMatrix_Private(DS ds,DSMatType dst,DSMatType src)
241: {
243: PetscInt j,m,off,ld;
244: PetscScalar *S,*D;
247: ld = ds->ld;
248: m = ds->n-ds->l;
249: off = ds->l+ds->l*ld;
250: S = ds->mat[src];
251: D = ds->mat[dst];
252: for (j=0;j<m;j++) {
253: PetscMemcpy(D+off+j*ld,S+off+j*ld,m*sizeof(PetscScalar));
254: }
255: return(0);
256: }
260: PetscErrorCode DSPermuteColumns_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat,PetscInt *perm)
261: {
262: PetscInt i,j,k,p,ld;
263: PetscScalar *Q,rtmp;
266: ld = ds->ld;
267: Q = ds->mat[mat];
268: for (i=l;i<n;i++) {
269: p = perm[i];
270: if (p != i) {
271: j = i + 1;
272: while (perm[j] != i) j++;
273: perm[j] = p; perm[i] = i;
274: /* swap columns i and j */
275: for (k=0;k<n;k++) {
276: rtmp = Q[k+p*ld]; Q[k+p*ld] = Q[k+i*ld]; Q[k+i*ld] = rtmp;
277: }
278: }
279: }
280: return(0);
281: }
285: PetscErrorCode DSPermuteRows_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat,PetscInt *perm)
286: {
287: PetscInt i,j,m=ds->m,k,p,ld;
288: PetscScalar *Q,rtmp;
291: if (m==0) SETERRQ(((PetscObject)ds)->comm,PETSC_ERR_ARG_WRONG,"m was not set");
292: ld = ds->ld;
293: Q = ds->mat[mat];
294: for (i=l;i<n;i++) {
295: p = perm[i];
296: if (p != i) {
297: j = i + 1;
298: while (perm[j] != i) j++;
299: perm[j] = p; perm[i] = i;
300: /* swap rows i and j */
301: for (k=0;k<m;k++) {
302: rtmp = Q[p+k*ld]; Q[p+k*ld] = Q[i+k*ld]; Q[i+k*ld] = rtmp;
303: }
304: }
305: }
306: return(0);
307: }
311: PetscErrorCode DSPermuteBoth_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat1,DSMatType mat2,PetscInt *perm)
312: {
313: PetscInt i,j,m=ds->m,k,p,ld;
314: PetscScalar *U,*VT,rtmp;
317: if (m==0) SETERRQ(((PetscObject)ds)->comm,PETSC_ERR_ARG_WRONG,"m was not set");
318: ld = ds->ld;
319: U = ds->mat[mat1];
320: VT = ds->mat[mat2];
321: for (i=l;i<n;i++) {
322: p = perm[i];
323: if (p != i) {
324: j = i + 1;
325: while (perm[j] != i) j++;
326: perm[j] = p; perm[i] = i;
327: /* swap columns i and j of U */
328: for (k=0;k<n;k++) {
329: rtmp = U[k+p*ld]; U[k+p*ld] = U[k+i*ld]; U[k+i*ld] = rtmp;
330: }
331: /* swap rows i and j of VT */
332: for (k=0;k<m;k++) {
333: rtmp = VT[p+k*ld]; VT[p+k*ld] = VT[i+k*ld]; VT[i+k*ld] = rtmp;
334: }
335: }
336: }
337: return(0);
338: }
342: /*
343: DSSetIdentity - Copy the identity (a diagonal matrix with ones) on the
344: active part of a matrix.
345: */
346: PetscErrorCode DSSetIdentity(DS ds,DSMatType mat)
347: {
349: PetscScalar *x;
350: PetscInt i,ld,n,l;
353: DSGetDimensions(ds,&n,PETSC_NULL,&l,PETSC_NULL);
354: DSGetLeadingDimension(ds,&ld);
355: DSGetArray(ds,mat,&x);
356: PetscLogEventBegin(DS_Other,ds,0,0,0);
357: PetscMemzero(&x[ld*l],ld*(n-l)*sizeof(PetscScalar));
358: for (i=l;i<n;i++) {
359: x[ld*i+i] = 1.0;
360: }
361: DSRestoreArray(ds,mat,&x);
362: PetscLogEventEnd(DS_Other,ds,0,0,0);
363: return(0);
364: }
368: /*
369: DSOrthogonalize - Orthogonalize the columns of a matrix.
371: Input Parameters:
372: + ds - the direct solver context
373: . mat - a matrix
374: - cols - number of columns to orthogonalize (starting from the column zero)
376: Output Parameter:
377: . lindcols - number of linearly independent columns of the matrix (can be PETSC_NULL)
378: */
379: PetscErrorCode DSOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscInt *lindcols)
380: {
381: #if defined(PETSC_MISSING_LAPACK_GEQRF) || defined(SLEPC_MISSING_LAPACK_ORGQR)
383: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQRF/ORGQR - Lapack routine is unavailable");
384: #else
385: PetscErrorCode ierr;
386: PetscInt n,l,ld;
387: PetscBLASInt ld_,rA,cA,info,ltau,lw;
388: PetscScalar *A,*tau,*w,saux;
394: DSGetDimensions(ds,&n,PETSC_NULL,&l,PETSC_NULL);
395: DSGetLeadingDimension(ds,&ld);
396: n = n - l;
397: if (cols > n) SETERRQ(((PetscObject)ds)->comm,PETSC_ERR_ARG_WRONG,"Invalid number of columns");
398: if (n == 0 || cols == 0) { return(0); }
399: DSGetArray(ds,mat,&A);
400: ltau = PetscBLASIntCast(PetscMin(cols,n));
401: ld_ = PetscBLASIntCast(ld);
402: rA = PetscBLASIntCast(n);
403: cA = PetscBLASIntCast(cols);
404: lw = -1;
405: LAPACKgeqrf_(&rA,&cA,A,&ld_,PETSC_NULL,&saux,&lw,&info);
406: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
407: lw = (PetscBLASInt)PetscRealPart(saux);
408: DSAllocateWork_Private(ds,lw+ltau,0,0);
409: tau = ds->work;
410: w = &tau[ltau];
411: PetscLogEventBegin(DS_Other,ds,0,0,0);
412: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
413: LAPACKgeqrf_(&rA,&cA,&A[ld*l+l],&ld_,tau,w,&lw,&info);
414: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
415: LAPACKorgqr_(&rA,<au,<au,&A[ld*l+l],&ld_,tau,w,&lw,&info);
416: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);
417: PetscFPTrapPop();
418: PetscLogEventEnd(DS_Other,ds,0,0,0);
419: DSRestoreArray(ds,mat,&A);
420: PetscObjectStateIncrease((PetscObject)ds);
421: if (lindcols) *lindcols = ltau;
422: return(0);
423: #endif
424: }
428: /*
429: DSPseudoOrthogonalize - Orthogonalize the columns of a matrix with Modified
430: Gram-Schmidt in an indefinite inner product space defined by a signature.
432: Input Parameters:
433: + ds - the direct solver context
434: . mat - the matrix
435: . cols - number of columns to orthogonalize (starting from the column zero)
436: - s - the signature that defines the inner product
438: Output Parameter:
439: + lindcols - linear independent columns of the matrix (can be PETSC_NULL)
440: - ns - the new norm of the vectors (can be PETSC_NULL)
441: */
442: PetscErrorCode DSPseudoOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscReal *s,PetscInt *lindcols,PetscReal *ns)
443: {
444: PetscErrorCode ierr;
445: PetscInt i,j,k,l,n,ld;
446: PetscBLASInt one=1,rA_;
447: PetscScalar alpha,*A,*A_,*m,*h,nr0;
448: PetscReal nr_o,nr,*ns_;
456: DSGetDimensions(ds,&n,PETSC_NULL,&l,PETSC_NULL);
457: DSGetLeadingDimension(ds,&ld);
458: n = n - l;
459: if (cols > n) SETERRQ(((PetscObject)ds)->comm,PETSC_ERR_ARG_WRONG,"Invalid number of columns");
460: if (n == 0 || cols == 0) { return(0); }
461: rA_ = PetscBLASIntCast(n);
462: DSGetArray(ds,mat,&A_);
463: A = &A_[ld*l+l];
464: DSAllocateWork_Private(ds,n+cols,ns?0:cols,0);
465: m = ds->work;
466: h = &m[n];
467: ns_ = ns ? ns : ds->rwork;
468: PetscLogEventBegin(DS_Other,ds,0,0,0);
469: for (i=0; i<cols; i++) {
470: /* m <- diag(s)*A[i] */
471: for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
472: /* nr_o <- mynorm(A[i]'*m), mynorm(x) = sign(x)*sqrt(|x|) */
473: SlepcDenseMatProd(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
474: nr = nr_o = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
475: for (j=0; j<3 && i>0; j++) {
476: /* h <- A[0:i-1]'*m */
477: SlepcDenseMatProd(h,i,0.0,1.0,A,ld,n,i,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
478: /* h <- diag(ns)*h */
479: for (k=0; k<i; k++) h[k] *= ns_[k];
480: /* A[i] <- A[i] - A[0:i-1]*h */
481: SlepcDenseMatProd(&A[ld*i],ld,1.0,-1.0,A,ld,n,i,PETSC_FALSE,h,i,i,1,PETSC_FALSE);
482: /* m <- diag(s)*A[i] */
483: for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
484: /* nr_o <- mynorm(A[i]'*m) */
485: SlepcDenseMatProd(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
486: nr = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
487: if (PetscAbs(nr) < PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_SELF,1, "Linear dependency detected");
488: if (PetscAbs(nr) > 0.7*PetscAbs(nr_o)) break;
489: nr_o = nr;
490: }
491: ns_[i] = PetscSign(nr);
492: /* A[i] <- A[i]/|nr| */
493: alpha = 1.0/PetscAbs(nr);
494: BLASscal_(&rA_,&alpha,&A[i*ld],&one);
495: }
496: PetscLogEventEnd(DS_Other,ds,0,0,0);
497: DSRestoreArray(ds,mat,&A_);
498: PetscObjectStateIncrease((PetscObject)ds);
499: if (lindcols) *lindcols = cols;
500: return(0);
501: }