1: /*
2: DS operations: DSSolve(), DSVectors(), etc
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2014, 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*/
28: /*@
29: DSGetLeadingDimension - Returns the leading dimension of the allocated
30: matrices.
32: Not Collective
34: Input Parameter:
35: . ds - the direct solver context
37: Output Parameter:
38: . ld - leading dimension (maximum allowed dimension for the matrices)
40: Level: advanced
42: .seealso: DSAllocate(), DSSetDimensions()
43: @*/
44: PetscErrorCode DSGetLeadingDimension(DS ds,PetscInt *ld) 45: {
49: *ld = ds->ld;
50: return(0);
51: }
55: /*@
56: DSSetState - Change the state of the DS object.
58: Logically Collective on DS 60: Input Parameters:
61: + ds - the direct solver context
62: - state - the new state
64: Notes:
65: The state indicates that the dense system is in an initial state (raw),
66: in an intermediate state (such as tridiagonal, Hessenberg or
67: Hessenberg-triangular), in a condensed state (such as diagonal, Schur or
68: generalized Schur), or in a truncated state.
70: This function is normally used to return to the raw state when the
71: condensed structure is destroyed.
73: Level: advanced
75: .seealso: DSGetState()
76: @*/
77: PetscErrorCode DSSetState(DS ds,DSStateType state) 78: {
84: switch (state) {
85: case DS_STATE_RAW:
86: case DS_STATE_INTERMEDIATE:
87: case DS_STATE_CONDENSED:
88: case DS_STATE_TRUNCATED:
89: if (ds->state<state) { PetscInfo(ds,"DS state has been increased\n"); }
90: ds->state = state;
91: break;
92: default: 93: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Wrong state");
94: }
95: return(0);
96: }
100: /*@
101: DSGetState - Returns the current state.
103: Not Collective
105: Input Parameter:
106: . ds - the direct solver context
108: Output Parameter:
109: . state - current state
111: Level: advanced
113: .seealso: DSSetState()
114: @*/
115: PetscErrorCode DSGetState(DS ds,DSStateType *state)116: {
120: *state = ds->state;
121: return(0);
122: }
126: /*@
127: DSSetDimensions - Resize the matrices in the DS object.
129: Logically Collective on DS131: Input Parameters:
132: + ds - the direct solver context
133: . n - the new size
134: . m - the new column size (only for DSSVD)
135: . l - number of locked (inactive) leading columns
136: - k - intermediate dimension (e.g., position of arrow)
138: Notes:
139: The internal arrays are not reallocated.
141: The value m is not used except in the case of DSSVD, pass 0 otherwise.
143: Level: intermediate
145: .seealso: DSGetDimensions(), DSAllocate()
146: @*/
147: PetscErrorCode DSSetDimensions(DS ds,PetscInt n,PetscInt m,PetscInt l,PetscInt k)148: {
151: DSCheckAlloc(ds,1);
156: if (n==PETSC_DECIDE || n==PETSC_DEFAULT) {
157: ds->n = ds->ld;
158: } else {
159: if (n<0 || n>ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of n. Must be between 0 and ld");
160: if (ds->extrarow && n+1>ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"A value of n equal to ld leaves no room for extra row");
161: ds->n = n;
162: }
163: ds->t = ds->n; /* truncated length equal to the new dimension */
164: if (m) {
165: if (m==PETSC_DECIDE || m==PETSC_DEFAULT) {
166: ds->m = ds->ld;
167: } else {
168: if (m<0 || m>ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of m. Must be between 0 and ld");
169: ds->m = m;
170: }
171: }
172: if (l==PETSC_DECIDE || l==PETSC_DEFAULT) {
173: ds->l = 0;
174: } else {
175: if (l<0 || l>ds->n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of l. Must be between 0 and n");
176: ds->l = l;
177: }
178: if (k==PETSC_DECIDE || k==PETSC_DEFAULT) {
179: ds->k = ds->n/2;
180: } else {
181: if (k<0 || k>ds->n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of k. Must be between 0 and n");
182: ds->k = k;
183: }
184: return(0);
185: }
189: /*@
190: DSGetDimensions - Returns the current dimensions.
192: Not Collective
194: Input Parameter:
195: . ds - the direct solver context
197: Output Parameter:
198: + n - the current size
199: . m - the current column size (only for DSSVD)
200: . l - number of locked (inactive) leading columns
201: . k - intermediate dimension (e.g., position of arrow)
202: - t - truncated length
204: Level: intermediate
206: Note:
207: The t parameter makes sense only if DSTruncate() has been called.
208: Otherwise its value equals n.
210: .seealso: DSSetDimensions(), DSTruncate()
211: @*/
212: PetscErrorCode DSGetDimensions(DS ds,PetscInt *n,PetscInt *m,PetscInt *l,PetscInt *k,PetscInt *t)213: {
216: DSCheckAlloc(ds,1);
217: if (n) *n = ds->n;
218: if (m) *m = ds->m;
219: if (l) *l = ds->l;
220: if (k) *k = ds->k;
221: if (t) *t = ds->t;
222: return(0);
223: }
227: /*@
228: DSTruncate - Truncates the system represented in the DS object.
230: Logically Collective on DS232: Input Parameters:
233: + ds - the direct solver context
234: - n - the new size
236: Note:
237: The new size is set to n. In cases where the extra row is meaningful,
238: the first n elements are kept as the extra row for the new system.
240: Level: advanced
242: .seealso: DSSetDimensions(), DSSetExtraRow(), DSStateType243: @*/
244: PetscErrorCode DSTruncate(DS ds,PetscInt n)245: {
250: DSCheckAlloc(ds,1);
251: DSCheckSolved(ds,1);
253: if (!ds->ops->truncate) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
254: if (n<ds->l || n>ds->n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of n. Must be between l and n");
255: PetscLogEventBegin(DS_Other,ds,0,0,0);
256: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
257: (*ds->ops->truncate)(ds,n);
258: PetscFPTrapPop();
259: PetscLogEventEnd(DS_Other,ds,0,0,0);
260: ds->state = DS_STATE_TRUNCATED;
261: PetscObjectStateIncrease((PetscObject)ds);
262: return(0);
263: }
267: /*@
268: DSGetMat - Returns a sequential dense Mat object containing the requested
269: matrix.
271: Not Collective
273: Input Parameters:
274: + ds - the direct solver context
275: - m - the requested matrix
277: Output Parameter:
278: . A - Mat object
280: Notes:
281: The Mat is created with sizes equal to the current DS dimensions (nxm),
282: then it is filled with the values that would be obtained with DSGetArray()
283: (not DSGetArrayReal()). If the DS was truncated, then the number of rows
284: is equal to the dimension prior to truncation, see DSTruncate().
285: The communicator is always PETSC_COMM_SELF.
287: When no longer needed, the user can either destroy the matrix or call
288: DSRestoreMat(). The latter will copy back the modified values.
290: Level: advanced
292: .seealso: DSRestoreMat(), DSSetDimensions(), DSGetArray(), DSGetArrayReal(), DSTruncate()
293: @*/
294: PetscErrorCode DSGetMat(DS ds,DSMatType m,Mat *A)295: {
297: PetscInt j,rows,cols,arows,acols;
298: PetscBool create=PETSC_FALSE;
299: PetscScalar *pA,*M;
303: DSCheckAlloc(ds,1);
306: if (m>=DS_NUM_MAT) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid matrix");
307: if (!ds->mat[m]) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Requested matrix was not created in this DS");
309: rows = PetscMax(ds->n,ds->t);
310: cols = ds->m? ds->m: ds->n;
311: if (!ds->omat[m]) create=PETSC_TRUE;
312: else {
313: MatGetSize(ds->omat[m],&arows,&acols);
314: if (arows!=rows || acols!=cols) {
315: MatDestroy(&ds->omat[m]);
316: create=PETSC_TRUE;
317: }
318: }
319: if (create) {
320: MatCreateSeqDense(PETSC_COMM_SELF,rows,cols,NULL,&ds->omat[m]);
321: }
322: PetscObjectReference((PetscObject)ds->omat[m]);
323: *A = ds->omat[m];
324: M = ds->mat[m];
325: MatDenseGetArray(*A,&pA);
326: for (j=0;j<cols;j++) {
327: PetscMemcpy(pA+j*rows,M+j*ds->ld,rows*sizeof(PetscScalar));
328: }
329: MatDenseRestoreArray(*A,&pA);
330: return(0);
331: }
335: /*@
336: DSRestoreMat - Restores the matrix after DSGetMat() was called.
338: Not Collective
340: Input Parameters:
341: + ds - the direct solver context
342: . m - the requested matrix
343: - A - the fetched Mat object
345: Notes:
346: A call to this function must match a previous call of DSGetMat().
347: The effect is that the contents of the Mat are copied back to the
348: DS internal array, and the matrix is destroyed.
350: It is not compulsory to call this function, the matrix obtained with
351: DSGetMat() can simply be destroyed if entries need not be copied back.
353: Level: advanced
355: .seealso: DSGetMat(), DSRestoreArray(), DSRestoreArrayReal()
356: @*/
357: PetscErrorCode DSRestoreMat(DS ds,DSMatType m,Mat *A)358: {
360: PetscInt j,rows,cols;
361: PetscScalar *pA,*M;
365: DSCheckAlloc(ds,1);
368: if (m>=DS_NUM_MAT) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid matrix");
369: if (!ds->omat[m]) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"DSRestoreMat must match a previous call to DSGetMat");
370: if (ds->omat[m]!=*A) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Mat argument is not the same as the one obtained with DSGetMat");
372: MatGetSize(*A,&rows,&cols);
373: M = ds->mat[m];
374: MatDenseGetArray(*A,&pA);
375: for (j=0;j<cols;j++) {
376: PetscMemcpy(M+j*ds->ld,pA+j*rows,rows*sizeof(PetscScalar));
377: }
378: MatDenseRestoreArray(*A,&pA);
379: MatDestroy(A);
380: return(0);
381: }
385: /*@C
386: DSGetArray - Returns a pointer to one of the internal arrays used to
387: represent matrices. You MUST call DSRestoreArray() when you no longer
388: need to access the array.
390: Not Collective
392: Input Parameters:
393: + ds - the direct solver context
394: - m - the requested matrix
396: Output Parameter:
397: . a - pointer to the values
399: Level: advanced
401: .seealso: DSRestoreArray(), DSGetArrayReal()
402: @*/
403: PetscErrorCode DSGetArray(DS ds,DSMatType m,PetscScalar *a[])404: {
407: DSCheckAlloc(ds,1);
409: if (m>=DS_NUM_MAT) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid matrix");
410: if (!ds->mat[m]) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Requested matrix was not created in this DS");
411: *a = ds->mat[m];
412: CHKMEMQ;
413: return(0);
414: }
418: /*@C
419: DSRestoreArray - Restores the matrix after DSGetArray() was called.
421: Not Collective
423: Input Parameters:
424: + ds - the direct solver context
425: . m - the requested matrix
426: - a - pointer to the values
428: Level: advanced
430: .seealso: DSGetArray(), DSGetArrayReal()
431: @*/
432: PetscErrorCode DSRestoreArray(DS ds,DSMatType m,PetscScalar *a[])433: {
438: DSCheckAlloc(ds,1);
440: if (m>=DS_NUM_MAT) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid matrix");
441: CHKMEMQ;
442: *a = 0;
443: PetscObjectStateIncrease((PetscObject)ds);
444: return(0);
445: }
449: /*@C
450: DSGetArrayReal - Returns a pointer to one of the internal arrays used to
451: represent real matrices. You MUST call DSRestoreArrayReal() when you no longer
452: need to access the array.
454: Not Collective
456: Input Parameters:
457: + ds - the direct solver context
458: - m - the requested matrix
460: Output Parameter:
461: . a - pointer to the values
463: Level: advanced
465: .seealso: DSRestoreArrayReal(), DSGetArray()
466: @*/
467: PetscErrorCode DSGetArrayReal(DS ds,DSMatType m,PetscReal *a[])468: {
471: DSCheckAlloc(ds,1);
473: if (m>=DS_NUM_MAT) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid matrix");
474: if (!ds->rmat[m]) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Requested matrix was not created in this DS");
475: *a = ds->rmat[m];
476: CHKMEMQ;
477: return(0);
478: }
482: /*@C
483: DSRestoreArrayReal - Restores the matrix after DSGetArrayReal() was called.
485: Not Collective
487: Input Parameters:
488: + ds - the direct solver context
489: . m - the requested matrix
490: - a - pointer to the values
492: Level: advanced
494: .seealso: DSGetArrayReal(), DSGetArray()
495: @*/
496: PetscErrorCode DSRestoreArrayReal(DS ds,DSMatType m,PetscReal *a[])497: {
502: DSCheckAlloc(ds,1);
504: if (m>=DS_NUM_MAT) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid matrix");
505: CHKMEMQ;
506: *a = 0;
507: PetscObjectStateIncrease((PetscObject)ds);
508: return(0);
509: }
513: /*@
514: DSSolve - Solves the problem.
516: Logically Collective on DS518: Input Parameters:
519: + ds - the direct solver context
520: . eigr - array to store the computed eigenvalues (real part)
521: - eigi - array to store the computed eigenvalues (imaginary part)
523: Note:
524: This call brings the dense system to condensed form. No ordering
525: of the eigenvalues is enforced (for this, call DSSort() afterwards).
527: Level: intermediate
529: .seealso: DSSort(), DSStateType530: @*/
531: PetscErrorCode DSSolve(DS ds,PetscScalar *eigr,PetscScalar *eigi)532: {
537: DSCheckAlloc(ds,1);
539: if (ds->state>=DS_STATE_CONDENSED) return(0);
540: if (!ds->ops->solve[ds->method]) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"The specified method number does not exist for this DS");
541: PetscLogEventBegin(DS_Solve,ds,0,0,0);
542: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
543: (*ds->ops->solve[ds->method])(ds,eigr,eigi);
544: PetscFPTrapPop();
545: PetscLogEventEnd(DS_Solve,ds,0,0,0);
546: ds->state = DS_STATE_CONDENSED;
547: PetscObjectStateIncrease((PetscObject)ds);
548: return(0);
549: }
553: /*@
554: DSComputeFunction - Computes a matrix function.
556: Logically Collective on DS558: Input Parameters:
559: + ds - the direct solver context
560: - f - the function to evaluate
562: Note:
563: This function evaluates F = f(A), where the input and the result matrices
564: are stored in DS_MAT_A and DS_MAT_F, respectively.
566: Level: intermediate
568: .seealso: DSSetFunctionMethod(), DSMatType, SlepcFunction
569: @*/
570: PetscErrorCode DSComputeFunction(DS ds,SlepcFunction f)571: {
576: DSCheckAlloc(ds,1);
578: if (!ds->ops->computefun[f][ds->funmethod]) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"The specified function method number does not exist for this DS");
579: if (!ds->mat[DS_MAT_F]) { DSAllocateMat_Private(ds,DS_MAT_F); }
580: PetscLogEventBegin(DS_Function,ds,0,0,0);
581: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
582: (*ds->ops->computefun[f][ds->funmethod])(ds);
583: PetscFPTrapPop();
584: PetscLogEventEnd(DS_Function,ds,0,0,0);
585: PetscObjectStateIncrease((PetscObject)ds);
586: return(0);
587: }
591: /*@
592: DSSort - Sorts the result of DSSolve() according to a given sorting
593: criterion.
595: Logically Collective on DS597: Input Parameters:
598: + ds - the direct solver context
599: . eigr - array containing the computed eigenvalues (real part)
600: . eigi - array containing the computed eigenvalues (imaginary part)
601: . rr - (optional) array containing auxiliary values (real part)
602: - ri - (optional) array containing auxiliary values (imaginary part)
604: Input/Output Parameter:
605: . k - (optional) number of elements in the leading group
607: Notes:
608: This routine sorts the arrays provided in eigr and eigi, and also
609: sorts the dense system stored inside ds (assumed to be in condensed form).
610: The sorting criterion is specified with DSSetSlepcSC().
612: If arrays rr and ri are provided, then a (partial) reordering based on these
613: values rather than on the eigenvalues is performed. In symmetric problems
614: a total order is obtained (parameter k is ignored), but otherwise the result
615: is sorted only partially. In this latter case, it is only guaranteed that
616: all the first k elements satisfy the comparison with any of the last n-k
617: elements. The output value of parameter k is the final number of elements in
618: the first set.
620: Level: intermediate
622: .seealso: DSSolve(), DSSetSlepcSC()
623: @*/
624: PetscErrorCode DSSort(DS ds,PetscScalar *eigr,PetscScalar *eigi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)625: {
630: DSCheckSolved(ds,1);
633: if (ds->state==DS_STATE_TRUNCATED) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Cannot sort a truncated DS");
634: if (!ds->ops->sort) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
635: if (!ds->sc) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must provide a sorting criterion first");
636: if (k && !rr) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Argument k can only be used together with rr");
637: PetscLogEventBegin(DS_Other,ds,0,0,0);
638: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
639: (*ds->ops->sort)(ds,eigr,eigi,rr,ri,k);
640: PetscFPTrapPop();
641: PetscLogEventEnd(DS_Other,ds,0,0,0);
642: PetscObjectStateIncrease((PetscObject)ds);
643: return(0);
644: }
648: /*@
649: DSVectors - Compute vectors associated to the dense system such
650: as eigenvectors.
652: Logically Collective on DS654: Input Parameters:
655: + ds - the direct solver context
656: - mat - the matrix, used to indicate which vectors are required
658: Input/Output Parameter:
659: - j - (optional) index of vector to be computed
661: Output Parameter:
662: . rnorm - (optional) computed residual norm
664: Notes:
665: Allowed values for mat are DS_MAT_X, DS_MAT_Y, DS_MAT_U and DS_MAT_VT, to
666: compute right or left eigenvectors, or left or right singular vectors,
667: respectively.
669: If NULL is passed in argument j then all vectors are computed,
670: otherwise j indicates which vector must be computed. In real non-symmetric
671: problems, on exit the index j will be incremented when a complex conjugate
672: pair is found.
674: This function can be invoked after the dense problem has been solved,
675: to get the residual norm estimate of the associated Ritz pair. In that
676: case, the relevant information is returned in rnorm.
678: For computing eigenvectors, LAPACK's _trevc is used so the matrix must
679: be in (quasi-)triangular form, or call DSSolve() first.
681: Level: intermediate
683: .seealso: DSSolve()
684: @*/
685: PetscErrorCode DSVectors(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)686: {
691: DSCheckAlloc(ds,1);
693: if (!ds->ops->vectors) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
694: if (rnorm && !j) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must give a value of j");
695: if (!ds->mat[mat]) { DSAllocateMat_Private(ds,mat); }
696: PetscLogEventBegin(DS_Vectors,ds,0,0,0);
697: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
698: (*ds->ops->vectors)(ds,mat,j,rnorm);
699: PetscFPTrapPop();
700: PetscLogEventEnd(DS_Vectors,ds,0,0,0);
701: PetscObjectStateIncrease((PetscObject)ds);
702: return(0);
703: }
707: /*@
708: DSNormalize - Normalize a column or all the columns of a matrix. Considers
709: the case when the columns represent the real and the imaginary part of a vector.
711: Logically Collective on DS713: Input Parameter:
714: + ds - the direct solver context
715: . mat - the matrix to be modified
716: - col - the column to normalize or -1 to normalize all of them
718: Notes:
719: The columns are normalized with respect to the 2-norm.
721: If col and col+1 (or col-1 and col) represent the real and the imaginary
722: part of a vector, both columns are scaled.
724: Level: advanced
726: .seealso: DSMatType727: @*/
728: PetscErrorCode DSNormalize(DS ds,DSMatType mat,PetscInt col)729: {
734: DSCheckSolved(ds,1);
737: if (!ds->ops->normalize) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
738: if (col<-1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"col should be at least minus one");
739: PetscLogEventBegin(DS_Other,ds,0,0,0);
740: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
741: (*ds->ops->normalize)(ds,mat,col);
742: PetscFPTrapPop();
743: PetscLogEventEnd(DS_Other,ds,0,0,0);
744: return(0);
745: }
749: /*@C
750: DSUpdateExtraRow - Performs all necessary operations so that the extra
751: row gets up-to-date after a call to DSSolve().
753: Not Collective
755: Input Parameters:
756: . ds - the direct solver context
758: Level: advanced
760: .seealso: DSSolve(), DSSetExtraRow()
761: @*/
762: PetscErrorCode DSUpdateExtraRow(DS ds)763: {
768: DSCheckAlloc(ds,1);
769: if (!ds->ops->update) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
770: if (!ds->extrarow) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Should have called DSSetExtraRow");
771: PetscLogEventBegin(DS_Other,ds,0,0,0);
772: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
773: (*ds->ops->update)(ds);
774: PetscFPTrapPop();
775: PetscLogEventEnd(DS_Other,ds,0,0,0);
776: return(0);
777: }
781: /*@C
782: DSCond - Compute the inf-norm condition number of the first matrix
783: as cond(A) = norm(A)*norm(inv(A)).
785: Not Collective
787: Input Parameters:
788: + ds - the direct solver context
789: - cond - the computed condition number
791: Level: advanced
793: .seealso: DSSolve()
794: @*/
795: PetscErrorCode DSCond(DS ds,PetscReal *cond)796: {
801: DSCheckAlloc(ds,1);
803: if (!ds->ops->cond) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
804: PetscLogEventBegin(DS_Other,ds,0,0,0);
805: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
806: (*ds->ops->cond)(ds,cond);
807: PetscFPTrapPop();
808: PetscLogEventEnd(DS_Other,ds,0,0,0);
809: return(0);
810: }
814: /*@C
815: DSTranslateHarmonic - Computes a translation of the dense system.
817: Logically Collective on DS819: Input Parameters:
820: + ds - the direct solver context
821: . tau - the translation amount
822: . beta - last component of vector b
823: - recover - boolean flag to indicate whether to recover or not
825: Output Parameters:
826: + g - the computed vector (optional)
827: - gamma - scale factor (optional)
829: Notes:
830: This function is intended for use in the context of Krylov methods only.
831: It computes a translation of a Krylov decomposition in order to extract
832: eigenpair approximations by harmonic Rayleigh-Ritz.
833: The matrix is updated as A + g*b' where g = (A-tau*eye(n))'\b and
834: vector b is assumed to be beta*e_n^T.
836: The gamma factor is defined as sqrt(1+g'*g) and can be interpreted as
837: the factor by which the residual of the Krylov decomposition is scaled.
839: If the recover flag is activated, the computed translation undoes the
840: translation done previously. In that case, parameter tau is ignored.
842: Level: developer
843: @*/
844: PetscErrorCode DSTranslateHarmonic(DS ds,PetscScalar tau,PetscReal beta,PetscBool recover,PetscScalar *g,PetscReal *gamma)845: {
850: DSCheckAlloc(ds,1);
851: if (!ds->ops->transharm) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
852: PetscLogEventBegin(DS_Other,ds,0,0,0);
853: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
854: (*ds->ops->transharm)(ds,tau,beta,recover,g,gamma);
855: PetscFPTrapPop();
856: PetscLogEventEnd(DS_Other,ds,0,0,0);
857: ds->state = DS_STATE_RAW;
858: PetscObjectStateIncrease((PetscObject)ds);
859: return(0);
860: }
864: /*@C
865: DSTranslateRKS - Computes a modification of the dense system corresponding
866: to an update of the shift in a rational Krylov method.
868: Logically Collective on DS870: Input Parameters:
871: + ds - the direct solver context
872: - alpha - the translation amount
874: Notes:
875: This function is intended for use in the context of Krylov methods only.
876: It takes the leading (k+1,k) submatrix of A, containing the truncated
877: Rayleigh quotient of a Krylov-Schur relation computed from a shift
878: sigma1 and transforms it to obtain a Krylov relation as if computed
879: from a different shift sigma2. The new matrix is computed as
880: 1.0/alpha*(eye(k)-Q*inv(R)), where [Q,R]=qr(eye(k)-alpha*A) and
881: alpha = sigma1-sigma2.
883: Matrix Q is placed in DS_MAT_Q so that it can be used to update the
884: Krylov basis.
886: Level: developer
887: @*/
888: PetscErrorCode DSTranslateRKS(DS ds,PetscScalar alpha)889: {
894: DSCheckAlloc(ds,1);
895: if (!ds->ops->transrks) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
896: PetscLogEventBegin(DS_Other,ds,0,0,0);
897: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
898: (*ds->ops->transrks)(ds,alpha);
899: PetscFPTrapPop();
900: PetscLogEventEnd(DS_Other,ds,0,0,0);
901: ds->state = DS_STATE_RAW;
902: ds->compact = PETSC_FALSE;
903: PetscObjectStateIncrease((PetscObject)ds);
904: return(0);
905: }