Actual source code: dspriv.c
1: /*
2: Private DS routines
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
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>
29: PetscErrorCode DSAllocateMat_Private(DS ds,DSMatType m)
30: {
31: PetscInt sz;
35: if (m==DS_MAT_T) sz = 3*ds->ld*sizeof(PetscScalar);
36: else if (m==DS_MAT_D) sz = ds->ld*sizeof(PetscScalar);
37: else sz = ds->ld*ds->ld*sizeof(PetscScalar);
38: if (ds->mat[m]) {
39: PetscFree(ds->mat[m]);
40: } else {
41: PetscLogObjectMemory(ds,sz);
42: }
43: PetscMalloc(sz,&ds->mat[m]);
44: PetscMemzero(ds->mat[m],sz);
45: return(0);
46: }
50: PetscErrorCode DSAllocateMatReal_Private(DS ds,DSMatType m)
51: {
52: PetscInt sz;
56: if (m==DS_MAT_T) sz = 3*ds->ld*sizeof(PetscReal);
57: else if (m==DS_MAT_D) sz = ds->ld*sizeof(PetscReal);
58: else sz = ds->ld*ds->ld*sizeof(PetscReal);
59: if (!ds->rmat[m]) {
60: PetscLogObjectMemory(ds,sz);
61: PetscMalloc(sz,&ds->rmat[m]);
62: }
63: PetscMemzero(ds->rmat[m],sz);
64: return(0);
65: }
69: PetscErrorCode DSAllocateWork_Private(DS ds,PetscInt s,PetscInt r,PetscInt i)
70: {
74: if (s>ds->lwork) {
75: PetscFree(ds->work);
76: PetscMalloc(s*sizeof(PetscScalar),&ds->work);
77: PetscLogObjectMemory(ds,(s-ds->lwork)*sizeof(PetscScalar));
78: ds->lwork = s;
79: }
80: if (r>ds->lrwork) {
81: PetscFree(ds->rwork);
82: PetscMalloc(r*sizeof(PetscReal),&ds->rwork);
83: PetscLogObjectMemory(ds,(r-ds->lrwork)*sizeof(PetscReal));
84: ds->lrwork = r;
85: }
86: if (i>ds->liwork) {
87: PetscFree(ds->iwork);
88: PetscMalloc(i*sizeof(PetscBLASInt),&ds->iwork);
89: PetscLogObjectMemory(ds,(i-ds->liwork)*sizeof(PetscBLASInt));
90: ds->liwork = i;
91: }
92: return(0);
93: }
97: PetscErrorCode DSViewMat_Private(DS ds,PetscViewer viewer,DSMatType m)
98: {
99: PetscErrorCode ierr;
100: PetscInt i,j,rows,cols;
101: PetscScalar *v;
102: PetscViewerFormat format;
103: #if defined(PETSC_USE_COMPLEX)
104: PetscBool allreal = PETSC_TRUE;
105: #endif
108: PetscViewerGetFormat(viewer,&format);
109: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
110: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
111: if (ds->state==DS_STATE_TRUNCATED && m>=DS_MAT_Q) rows = ds->t;
112: else rows = (m==DS_MAT_A && ds->extrarow)? ds->n+1: ds->n;
113: cols = (ds->m!=0)? ds->m: ds->n;
114: #if defined(PETSC_USE_COMPLEX)
115: /* determine if matrix has all real values */
116: v = ds->mat[m];
117: for (i=0;i<rows;i++)
118: for (j=0;j<cols;j++)
119: if (PetscImaginaryPart(v[i+j*ds->ld])) { allreal = PETSC_FALSE; break; }
120: #endif
121: if (format == PETSC_VIEWER_ASCII_MATLAB) {
122: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",rows,cols);
123: PetscViewerASCIIPrintf(viewer,"%s = [\n",DSMatName[m]);
124: } else {
125: PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[m]);
126: }
128: for (i=0;i<rows;i++) {
129: v = ds->mat[m]+i;
130: for (j=0;j<cols;j++) {
131: #if defined(PETSC_USE_COMPLEX)
132: if (allreal) {
133: PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));
134: } else {
135: PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei ",PetscRealPart(*v),PetscImaginaryPart(*v));
136: }
137: #else
138: PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);
139: #endif
140: v += ds->ld;
141: }
142: PetscViewerASCIIPrintf(viewer,"\n");
143: }
145: if (format == PETSC_VIEWER_ASCII_MATLAB) {
146: PetscViewerASCIIPrintf(viewer,"];\n");
147: }
148: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
149: PetscViewerFlush(viewer);
150: return(0);
151: }
155: PetscErrorCode DSSortEigenvalues_Private(DS ds,PetscScalar *wr,PetscScalar *wi,PetscInt *perm,PetscBool isghiep)
156: {
158: PetscScalar re,im,wi0;
159: PetscInt n,i,j,result,tmp1,tmp2=0,d=1;
162: n = ds->t; /* sort only first t pairs if truncated */
163: for (i=0;i<ds->n;i++) perm[i] = i;
164: /* insertion sort */
165: i=ds->l+1;
166: #if !defined(PETSC_USE_COMPLEX)
167: if (wi && wi[perm[i-1]]!=0.0) i++; /* initial value is complex */
168: #else
169: if (isghiep && PetscImaginaryPart(wr[perm[i-1]])!=0.0) i++;
170: #endif
171: for (;i<n;i+=d) {
172: re = wr[perm[i]];
173: if (wi) im = wi[perm[i]];
174: else im = 0.0;
175: tmp1 = perm[i];
176: #if !defined(PETSC_USE_COMPLEX)
177: if (im!=0.0) { d = 2; tmp2 = perm[i+1]; }
178: else d = 1;
179: #else
180: if (isghiep && PetscImaginaryPart(re)!=0.0) { d = 2; tmp2 = perm[i+1]; }
181: else d = 1;
182: #endif
183: j = i-1;
184: if (wi) wi0 = wi[perm[j]];
185: else wi0 = 0.0;
186: (*ds->comparison)(re,im,wr[perm[j]],wi0,&result,ds->comparisonctx);
187: while (result<0 && j>=ds->l) {
188: perm[j+d] = perm[j];
189: j--;
190: #if !defined(PETSC_USE_COMPLEX)
191: if (wi && wi[perm[j+1]]!=0)
192: #else
193: if (isghiep && PetscImaginaryPart(wr[perm[j+1]])!=0)
194: #endif
195: { perm[j+d] = perm[j]; j--; }
197: if (j>=ds->l) {
198: if (wi) wi0 = wi[perm[j]];
199: else wi0 = 0.0;
200: (*ds->comparison)(re,im,wr[perm[j]],wi0,&result,ds->comparisonctx);
201: }
202: }
203: perm[j+1] = tmp1;
204: if (d==2) perm[j+2] = tmp2;
205: }
206: return(0);
207: }
211: PetscErrorCode DSSortEigenvaluesReal_Private(DS ds,PetscReal *eig,PetscInt *perm)
212: {
214: PetscScalar re;
215: PetscInt i,j,result,tmp,l,n;
218: n = ds->t; /* sort only first t pairs if truncated */
219: l = ds->l;
220: for (i=0;i<n;i++) perm[i] = i;
221: /* insertion sort */
222: for (i=l+1;i<n;i++) {
223: re = eig[perm[i]];
224: j = i-1;
225: (*ds->comparison)(re,0.0,eig[perm[j]],0.0,&result,ds->comparisonctx);
226: while (result<0 && j>=l) {
227: tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
228: if (j>=l) {
229: (*ds->comparison)(re,0.0,eig[perm[j]],0.0,&result,ds->comparisonctx);
230: }
231: }
232: }
233: return(0);
234: }
238: /*
239: DSCopyMatrix_Private - Copies the trailing block of a matrix (from
240: rows/columns l to n).
241: */
242: PetscErrorCode DSCopyMatrix_Private(DS ds,DSMatType dst,DSMatType src)
243: {
245: PetscInt j,m,off,ld;
246: PetscScalar *S,*D;
249: ld = ds->ld;
250: m = ds->n-ds->l;
251: off = ds->l+ds->l*ld;
252: S = ds->mat[src];
253: D = ds->mat[dst];
254: for (j=0;j<m;j++) {
255: PetscMemcpy(D+off+j*ld,S+off+j*ld,m*sizeof(PetscScalar));
256: }
257: return(0);
258: }
262: PetscErrorCode DSPermuteColumns_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat,PetscInt *perm)
263: {
264: PetscInt i,j,k,p,ld;
265: PetscScalar *Q,rtmp;
268: ld = ds->ld;
269: Q = ds->mat[mat];
270: for (i=l;i<n;i++) {
271: p = perm[i];
272: if (p != i) {
273: j = i + 1;
274: while (perm[j] != i) j++;
275: perm[j] = p; perm[i] = i;
276: /* swap columns i and j */
277: for (k=0;k<n;k++) {
278: rtmp = Q[k+p*ld]; Q[k+p*ld] = Q[k+i*ld]; Q[k+i*ld] = rtmp;
279: }
280: }
281: }
282: return(0);
283: }
287: PetscErrorCode DSPermuteRows_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat,PetscInt *perm)
288: {
289: PetscInt i,j,m=ds->m,k,p,ld;
290: PetscScalar *Q,rtmp;
293: if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
294: ld = ds->ld;
295: Q = ds->mat[mat];
296: for (i=l;i<n;i++) {
297: p = perm[i];
298: if (p != i) {
299: j = i + 1;
300: while (perm[j] != i) j++;
301: perm[j] = p; perm[i] = i;
302: /* swap rows i and j */
303: for (k=0;k<m;k++) {
304: rtmp = Q[p+k*ld]; Q[p+k*ld] = Q[i+k*ld]; Q[i+k*ld] = rtmp;
305: }
306: }
307: }
308: return(0);
309: }
313: PetscErrorCode DSPermuteBoth_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat1,DSMatType mat2,PetscInt *perm)
314: {
315: PetscInt i,j,m=ds->m,k,p,ld;
316: PetscScalar *U,*VT,rtmp;
319: if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
320: ld = ds->ld;
321: U = ds->mat[mat1];
322: VT = ds->mat[mat2];
323: for (i=l;i<n;i++) {
324: p = perm[i];
325: if (p != i) {
326: j = i + 1;
327: while (perm[j] != i) j++;
328: perm[j] = p; perm[i] = i;
329: /* swap columns i and j of U */
330: for (k=0;k<n;k++) {
331: rtmp = U[k+p*ld]; U[k+p*ld] = U[k+i*ld]; U[k+i*ld] = rtmp;
332: }
333: /* swap rows i and j of VT */
334: for (k=0;k<m;k++) {
335: rtmp = VT[p+k*ld]; VT[p+k*ld] = VT[i+k*ld]; VT[i+k*ld] = rtmp;
336: }
337: }
338: }
339: return(0);
340: }
344: /*
345: DSSetIdentity - Copy the identity (a diagonal matrix with ones) on the
346: active part of a matrix.
347: */
348: PetscErrorCode DSSetIdentity(DS ds,DSMatType mat)
349: {
351: PetscScalar *x;
352: PetscInt i,ld,n,l;
355: DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
356: DSGetLeadingDimension(ds,&ld);
357: PetscLogEventBegin(DS_Other,ds,0,0,0);
358: DSGetArray(ds,mat,&x);
359: PetscMemzero(&x[ld*l],ld*(n-l)*sizeof(PetscScalar));
360: for (i=l;i<n;i++) {
361: x[ld*i+i] = 1.0;
362: }
363: DSRestoreArray(ds,mat,&x);
364: PetscLogEventEnd(DS_Other,ds,0,0,0);
365: return(0);
366: }
370: /*
371: DSComputeMatrix - Build the matrix associated to a nonlinear operator
372: T(lambda) or its derivative T'(lambda), given the parameter lambda, where
373: T(lambda) = sum_i E_i*f_i(lambda). The result is written in mat.
374: */
375: PetscErrorCode DSComputeMatrix(DS ds,PetscScalar lambda,PetscBool deriv,DSMatType mat)
376: {
378: PetscScalar *T,*E,alpha;
379: PetscInt i,ld,n;
380: PetscBLASInt k,inc=1;
383: DSGetDimensions(ds,&n,NULL,NULL,NULL,NULL);
384: DSGetLeadingDimension(ds,&ld);
385: PetscBLASIntCast(ld*n,&k);
386: PetscLogEventBegin(DS_Other,ds,0,0,0);
387: DSGetArray(ds,mat,&T);
388: PetscMemzero(T,k*sizeof(PetscScalar));
389: for (i=0;i<ds->nf;i++) {
390: if (deriv) {
391: FNEvaluateDerivative(ds->f[i],lambda,&alpha);
392: } else {
393: FNEvaluateFunction(ds->f[i],lambda,&alpha);
394: }
395: E = ds->mat[DSMatExtra[i]];
396: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k,&alpha,E,&inc,T,&inc));
397: }
398: DSRestoreArray(ds,mat,&T);
399: PetscLogEventEnd(DS_Other,ds,0,0,0);
400: return(0);
401: }
405: /*
406: DSOrthogonalize - Orthogonalize the columns of a matrix.
408: Input Parameters:
409: + ds - the direct solver context
410: . mat - a matrix
411: - cols - number of columns to orthogonalize (starting from the column zero)
413: Output Parameter:
414: . lindcols - number of linearly independent columns of the matrix (can be NULL)
415: */
416: PetscErrorCode DSOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscInt *lindcols)
417: {
418: #if defined(PETSC_MISSING_LAPACK_GEQRF) || defined(SLEPC_MISSING_LAPACK_ORGQR)
420: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQRF/ORGQR - Lapack routine is unavailable");
421: #else
422: PetscErrorCode ierr;
423: PetscInt n,l,ld;
424: PetscBLASInt ld_,rA,cA,info,ltau,lw;
425: PetscScalar *A,*tau,*w,saux;
431: DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
432: DSGetLeadingDimension(ds,&ld);
433: n = n - l;
434: if (cols > n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns");
435: if (n == 0 || cols == 0) return(0);
436: DSGetArray(ds,mat,&A);
437: PetscBLASIntCast(PetscMin(cols,n),<au);
438: PetscBLASIntCast(ld,&ld_);
439: PetscBLASIntCast(n,&rA);
440: PetscBLASIntCast(cols,&cA);
441: lw = -1;
442: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,A,&ld_,NULL,&saux,&lw,&info));
443: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
444: lw = (PetscBLASInt)PetscRealPart(saux);
445: DSAllocateWork_Private(ds,lw+ltau,0,0);
446: tau = ds->work;
447: w = &tau[ltau];
448: PetscLogEventBegin(DS_Other,ds,0,0,0);
449: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
450: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,&A[ld*l+l],&ld_,tau,w,&lw,&info));
451: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
452: PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&rA,<au,<au,&A[ld*l+l],&ld_,tau,w,&lw,&info));
453: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);
454: PetscFPTrapPop();
455: PetscLogEventEnd(DS_Other,ds,0,0,0);
456: DSRestoreArray(ds,mat,&A);
457: PetscObjectStateIncrease((PetscObject)ds);
458: if (lindcols) *lindcols = ltau;
459: return(0);
460: #endif
461: }
465: /*
466: DSPseudoOrthogonalize - Orthogonalize the columns of a matrix with Modified
467: Gram-Schmidt in an indefinite inner product space defined by a signature.
469: Input Parameters:
470: + ds - the direct solver context
471: . mat - the matrix
472: . cols - number of columns to orthogonalize (starting from the column zero)
473: - s - the signature that defines the inner product
475: Output Parameter:
476: + lindcols - linear independent columns of the matrix (can be NULL)
477: - ns - the new norm of the vectors (can be NULL)
478: */
479: PetscErrorCode DSPseudoOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscReal *s,PetscInt *lindcols,PetscReal *ns)
480: {
481: PetscErrorCode ierr;
482: PetscInt i,j,k,l,n,ld;
483: PetscBLASInt one=1,rA_;
484: PetscScalar alpha,*A,*A_,*m,*h,nr0;
485: PetscReal nr_o,nr,*ns_;
493: DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
494: DSGetLeadingDimension(ds,&ld);
495: n = n - l;
496: if (cols > n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns");
497: if (n == 0 || cols == 0) return(0);
498: PetscBLASIntCast(n,&rA_);
499: DSGetArray(ds,mat,&A_);
500: A = &A_[ld*l+l];
501: DSAllocateWork_Private(ds,n+cols,ns?0:cols,0);
502: m = ds->work;
503: h = &m[n];
504: ns_ = ns ? ns : ds->rwork;
505: PetscLogEventBegin(DS_Other,ds,0,0,0);
506: for (i=0; i<cols; i++) {
507: /* m <- diag(s)*A[i] */
508: for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
509: /* nr_o <- mynorm(A[i]'*m), mynorm(x) = sign(x)*sqrt(|x|) */
510: SlepcDenseMatProd(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
511: nr = nr_o = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
512: for (j=0; j<3 && i>0; j++) {
513: /* h <- A[0:i-1]'*m */
514: SlepcDenseMatProd(h,i,0.0,1.0,A,ld,n,i,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
515: /* h <- diag(ns)*h */
516: for (k=0; k<i; k++) h[k] *= ns_[k];
517: /* A[i] <- A[i] - A[0:i-1]*h */
518: SlepcDenseMatProd(&A[ld*i],ld,1.0,-1.0,A,ld,n,i,PETSC_FALSE,h,i,i,1,PETSC_FALSE);
519: /* m <- diag(s)*A[i] */
520: for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
521: /* nr_o <- mynorm(A[i]'*m) */
522: SlepcDenseMatProd(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
523: nr = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
524: if (PetscAbs(nr) < PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_SELF,1,"Linear dependency detected");
525: if (PetscAbs(nr) > 0.7*PetscAbs(nr_o)) break;
526: nr_o = nr;
527: }
528: ns_[i] = PetscSign(nr);
529: /* A[i] <- A[i]/|nr| */
530: alpha = 1.0/PetscAbs(nr);
531: PetscStackCallBLAS("BLASscal",BLASscal_(&rA_,&alpha,&A[i*ld],&one));
532: }
533: PetscLogEventEnd(DS_Other,ds,0,0,0);
534: DSRestoreArray(ds,mat,&A_);
535: PetscObjectStateIncrease((PetscObject)ds);
536: if (lindcols) *lindcols = cols;
537: return(0);
538: }