Actual source code: dsops.c
1: /*
2: DS operations: DSSolve(), DSVectors(), etc
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*/
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: {
48: if (ld) *ld = ds->ld;
49: return(0);
50: }
54: /*@
55: DSSetState - Change the state of the DS object.
57: Logically Collective on DS
59: Input Parameters:
60: + ds - the direct solver context
61: - state - the new state
63: Notes:
64: The state indicates that the dense system is in an initial state (raw),
65: in an intermediate state (such as tridiagonal, Hessenberg or
66: Hessenberg-triangular), in a condensed state (such as diagonal, Schur or
67: generalized Schur), or in a truncated state.
69: This function is normally used to return to the raw state when the
70: condensed structure is destroyed.
72: Level: advanced
74: .seealso: DSGetState()
75: @*/
76: PetscErrorCode DSSetState(DS ds,DSStateType state)
77: {
83: switch (state) {
84: case DS_STATE_RAW:
85: case DS_STATE_INTERMEDIATE:
86: case DS_STATE_CONDENSED:
87: case DS_STATE_TRUNCATED:
88: if (ds->state<state) { PetscInfo(ds,"DS state has been increased\n"); }
89: ds->state = state;
90: break;
91: default:
92: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Wrong state");
93: }
94: return(0);
95: }
99: /*@
100: DSGetState - Returns the current state.
102: Not Collective
104: Input Parameter:
105: . ds - the direct solver context
107: Output Parameter:
108: . state - current state
110: Level: advanced
112: .seealso: DSSetState()
113: @*/
114: PetscErrorCode DSGetState(DS ds,DSStateType *state)
115: {
118: if (state) *state = ds->state;
119: return(0);
120: }
124: /*@
125: DSSetDimensions - Resize the matrices in the DS object.
127: Logically Collective on DS
129: Input Parameters:
130: + ds - the direct solver context
131: . n - the new size
132: . m - the new column size (only for DSSVD)
133: . l - number of locked (inactive) leading columns
134: - k - intermediate dimension (e.g., position of arrow)
136: Notes:
137: The internal arrays are not reallocated.
139: The value m is not used except in the case of DSSVD, pass 0 otherwise.
141: Level: intermediate
143: .seealso: DSGetDimensions(), DSAllocate()
144: @*/
145: PetscErrorCode DSSetDimensions(DS ds,PetscInt n,PetscInt m,PetscInt l,PetscInt k)
146: {
153: if (!ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must call DSAllocate() first");
154: if (n) {
155: if (n==PETSC_DECIDE || n==PETSC_DEFAULT) {
156: ds->n = ds->ld;
157: } else {
158: if (n<1 || n>ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of n. Must be between 1 and ld");
159: 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");
160: ds->n = n;
161: }
162: ds->t = ds->n; /* truncated length equal to the new dimension */
163: }
164: if (m) {
165: if (m==PETSC_DECIDE || m==PETSC_DEFAULT) {
166: ds->m = ds->ld;
167: } else {
168: if (m<1 || m>ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of m. Must be between 1 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: if (n) *n = ds->n;
217: if (m) *m = ds->m;
218: if (l) *l = ds->l;
219: if (k) *k = ds->k;
220: if (t) *t = ds->t;
221: return(0);
222: }
226: /*@
227: DSTruncate - Truncates the system represented in the DS object.
229: Logically Collective on DS
231: Input Parameters:
232: + ds - the direct solver context
233: - n - the new size
235: Note:
236: The new size is set to n. In cases where the extra row is meaningful,
237: the first n elements are kept as the extra row for the new system.
239: Level: advanced
241: .seealso: DSSetDimensions(), DSSetExtraRow(), DSStateType
242: @*/
243: PetscErrorCode DSTruncate(DS ds,PetscInt n)
244: {
250: if (!ds->ops->truncate) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
251: if (ds->state<DS_STATE_CONDENSED) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must call DSSolve() first");
252: 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");
253: PetscLogEventBegin(DS_Other,ds,0,0,0);
254: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
255: (*ds->ops->truncate)(ds,n);
256: PetscFPTrapPop();
257: PetscLogEventEnd(DS_Other,ds,0,0,0);
258: ds->state = DS_STATE_TRUNCATED;
259: PetscObjectStateIncrease((PetscObject)ds);
260: return(0);
261: }
265: /*@C
266: DSGetArray - Returns a pointer to one of the internal arrays used to
267: represent matrices. You MUST call DSRestoreArray() when you no longer
268: need to access the array.
270: Not Collective
272: Input Parameters:
273: + ds - the direct solver context
274: - m - the requested matrix
276: Output Parameter:
277: . a - pointer to the values
279: Level: advanced
281: .seealso: DSRestoreArray(), DSGetArrayReal()
282: @*/
283: PetscErrorCode DSGetArray(DS ds,DSMatType m,PetscScalar *a[])
284: {
288: if (m<0 || m>=DS_NUM_MAT) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid matrix");
289: if (!ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must call DSAllocate() first");
290: if (!ds->mat[m]) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Requested matrix was not created in this DS");
291: *a = ds->mat[m];
292: CHKMEMQ;
293: return(0);
294: }
298: /*@C
299: DSRestoreArray - Restores the matrix after DSGetArray() was called.
301: Not Collective
303: Input Parameters:
304: + ds - the direct solver context
305: . m - the requested matrix
306: - a - pointer to the values
308: Level: advanced
310: .seealso: DSGetArray(), DSGetArrayReal()
311: @*/
312: PetscErrorCode DSRestoreArray(DS ds,DSMatType m,PetscScalar *a[])
313: {
319: if (m<0 || m>=DS_NUM_MAT) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid matrix");
320: CHKMEMQ;
321: *a = 0;
322: PetscObjectStateIncrease((PetscObject)ds);
323: return(0);
324: }
328: /*@C
329: DSGetArrayReal - Returns a pointer to one of the internal arrays used to
330: represent real matrices. You MUST call DSRestoreArray() when you no longer
331: need to access the array.
333: Not Collective
335: Input Parameters:
336: + ds - the direct solver context
337: - m - the requested matrix
339: Output Parameter:
340: . a - pointer to the values
342: Level: advanced
344: .seealso: DSRestoreArrayReal(), DSGetArray()
345: @*/
346: PetscErrorCode DSGetArrayReal(DS ds,DSMatType m,PetscReal *a[])
347: {
351: if (m<0 || m>=DS_NUM_MAT) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid matrix");
352: if (!ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must call DSAllocate() first");
353: if (!ds->rmat[m]) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Requested matrix was not created in this DS");
354: *a = ds->rmat[m];
355: CHKMEMQ;
356: return(0);
357: }
361: /*@C
362: DSRestoreArrayReal - Restores the matrix after DSGetArrayReal() was called.
364: Not Collective
366: Input Parameters:
367: + ds - the direct solver context
368: . m - the requested matrix
369: - a - pointer to the values
371: Level: advanced
373: .seealso: DSGetArrayReal(), DSGetArray()
374: @*/
375: PetscErrorCode DSRestoreArrayReal(DS ds,DSMatType m,PetscReal *a[])
376: {
382: if (m<0 || m>=DS_NUM_MAT) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid matrix");
383: CHKMEMQ;
384: *a = 0;
385: PetscObjectStateIncrease((PetscObject)ds);
386: return(0);
387: }
391: /*@
392: DSSolve - Solves the problem.
394: Logically Collective on DS
396: Input Parameters:
397: + ds - the direct solver context
398: . eigr - array to store the computed eigenvalues (real part)
399: - eigi - array to store the computed eigenvalues (imaginary part)
401: Note:
402: This call brings the dense system to condensed form. No ordering
403: of the eigenvalues is enforced (for this, call DSSort() afterwards).
405: Level: intermediate
407: .seealso: DSSort(), DSStateType
408: @*/
409: PetscErrorCode DSSolve(DS ds,PetscScalar *eigr,PetscScalar *eigi)
410: {
416: if (!ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must call DSAllocate() first");
417: if (ds->state>=DS_STATE_CONDENSED) return(0);
418: if (!ds->ops->solve[ds->method]) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"The specified method number does not exist for this DS");
419: PetscLogEventBegin(DS_Solve,ds,0,0,0);
420: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
421: (*ds->ops->solve[ds->method])(ds,eigr,eigi);
422: PetscFPTrapPop();
423: PetscLogEventEnd(DS_Solve,ds,0,0,0);
424: ds->state = DS_STATE_CONDENSED;
425: PetscObjectStateIncrease((PetscObject)ds);
426: return(0);
427: }
431: /*@
432: DSComputeFunction - Computes a matrix function.
434: Logically Collective on DS
436: Input Parameters:
437: + ds - the direct solver context
438: - f - the function to evaluate
440: Note:
441: This function evaluates F = f(A), where the input and the result matrices
442: are stored in DS_MAT_A and DS_MAT_F, respectively.
444: Level: intermediate
446: .seealso: DSSetFunctionMethod(), DSMatType, SlepcFunction
447: @*/
448: PetscErrorCode DSComputeFunction(DS ds,SlepcFunction f)
449: {
455: if (!ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must call DSAllocate() first");
456: 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");
457: if (!ds->mat[DS_MAT_F]) { DSAllocateMat_Private(ds,DS_MAT_F); }
458: PetscLogEventBegin(DS_Function,ds,0,0,0);
459: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
460: (*ds->ops->computefun[f][ds->funmethod])(ds);
461: PetscFPTrapPop();
462: PetscLogEventEnd(DS_Function,ds,0,0,0);
463: PetscObjectStateIncrease((PetscObject)ds);
464: return(0);
465: }
469: /*@
470: DSSort - Sorts the result of DSSolve() according to a given sorting
471: criterion.
473: Logically Collective on DS
475: Input Parameters:
476: + ds - the direct solver context
477: . eigr - array containing the computed eigenvalues (real part)
478: . eigi - array containing the computed eigenvalues (imaginary part)
479: . rr - (optional) array containing auxiliary values (real part)
480: - ri - (optional) array containing auxiliary values (imaginary part)
482: Input/Output Parameter:
483: . k - (optional) number of elements in the leading group
485: Notes:
486: This routine sorts the arrays provided in eigr and eigi, and also
487: sorts the dense system stored inside ds (assumed to be in condensed form).
488: The sorting criterion is specified with DSSetEigenvalueComparison().
490: If arrays rr and ri are provided, then a (partial) reordering based on these
491: values rather than on the eigenvalues is performed. In symmetric problems
492: a total order is obtained (parameter k is ignored), but otherwise the result
493: is sorted only partially. In this latter case, it is only guaranteed that
494: all the first k elements satisfy the comparison with any of the last n-k
495: elements. The output value of parameter k is the final number of elements in
496: the first set.
498: Level: intermediate
500: .seealso: DSSolve(), DSSetEigenvalueComparison()
501: @*/
502: PetscErrorCode DSSort(DS ds,PetscScalar *eigr,PetscScalar *eigi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
503: {
510: if (ds->state<DS_STATE_CONDENSED) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must call DSSolve() first");
511: if (ds->state==DS_STATE_TRUNCATED) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Cannot sort a truncated DS");
512: if (!ds->ops->sort) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
513: if (!ds->comparison) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must provide a sorting criterion with DSSetEigenvalueComparison() first");
514: if (k && !rr) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Argument k can only be used together with rr");
515: PetscLogEventBegin(DS_Other,ds,0,0,0);
516: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
517: (*ds->ops->sort)(ds,eigr,eigi,rr,ri,k);
518: PetscFPTrapPop();
519: PetscLogEventEnd(DS_Other,ds,0,0,0);
520: PetscObjectStateIncrease((PetscObject)ds);
521: return(0);
522: }
526: /*@
527: DSVectors - Compute vectors associated to the dense system such
528: as eigenvectors.
530: Logically Collective on DS
532: Input Parameters:
533: + ds - the direct solver context
534: - mat - the matrix, used to indicate which vectors are required
536: Input/Output Parameter:
537: - j - (optional) index of vector to be computed
539: Output Parameter:
540: . rnorm - (optional) computed residual norm
542: Notes:
543: Allowed values for mat are DS_MAT_X, DS_MAT_Y, DS_MAT_U and DS_MAT_VT, to
544: compute right or left eigenvectors, or left or right singular vectors,
545: respectively.
547: If NULL is passed in argument j then all vectors are computed,
548: otherwise j indicates which vector must be computed. In real non-symmetric
549: problems, on exit the index j will be incremented when a complex conjugate
550: pair is found.
552: This function can be invoked after the dense problem has been solved,
553: to get the residual norm estimate of the associated Ritz pair. In that
554: case, the relevant information is returned in rnorm.
556: For computing eigenvectors, LAPACK's _trevc is used so the matrix must
557: be in (quasi-)triangular form, or call DSSolve() first.
559: Level: intermediate
561: .seealso: DSSolve()
562: @*/
563: PetscErrorCode DSVectors(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
564: {
570: if (!ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must call DSAllocate() first");
571: if (!ds->ops->vectors) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
572: if (rnorm && !j) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must give a value of j");
573: if (!ds->mat[mat]) { DSAllocateMat_Private(ds,mat); }
574: PetscLogEventBegin(DS_Vectors,ds,0,0,0);
575: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
576: (*ds->ops->vectors)(ds,mat,j,rnorm);
577: PetscFPTrapPop();
578: PetscLogEventEnd(DS_Vectors,ds,0,0,0);
579: PetscObjectStateIncrease((PetscObject)ds);
580: return(0);
581: }
585: /*@
586: DSNormalize - Normalize a column or all the columns of a matrix. Considers
587: the case when the columns represent the real and the imaginary part of a vector.
589: Logically Collective on DS
591: Input Parameter:
592: + ds - the direct solver context
593: . mat - the matrix to be modified
594: - col - the column to normalize or -1 to normalize all of them
596: Notes:
597: The columns are normalized with respect to the 2-norm.
599: If col and col+1 (or col-1 and col) represent the real and the imaginary
600: part of a vector, both columns are scaled.
602: Level: advanced
604: .seealso: DSMatType
605: @*/
606: PetscErrorCode DSNormalize(DS ds,DSMatType mat,PetscInt col)
607: {
614: if (ds->state<DS_STATE_CONDENSED) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must call DSSolve() first");
615: if (!ds->ops->normalize) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
616: if (col<-1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"col should be at least minus one");
617: PetscLogEventBegin(DS_Other,ds,0,0,0);
618: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
619: (*ds->ops->normalize)(ds,mat,col);
620: PetscFPTrapPop();
621: PetscLogEventEnd(DS_Other,ds,0,0,0);
622: return(0);
623: }
627: /*@C
628: DSUpdateExtraRow - Performs all necessary operations so that the extra
629: row gets up-to-date after a call to DSSolve().
631: Not Collective
633: Input Parameters:
634: . ds - the direct solver context
636: Level: advanced
638: .seealso: DSSolve(), DSSetExtraRow()
639: @*/
640: PetscErrorCode DSUpdateExtraRow(DS ds)
641: {
646: if (!ds->ops->update) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
647: if (!ds->extrarow) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Should have called DSSetExtraRow");
648: if (!ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must call DSAllocate() first");
649: PetscLogEventBegin(DS_Other,ds,0,0,0);
650: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
651: (*ds->ops->update)(ds);
652: PetscFPTrapPop();
653: PetscLogEventEnd(DS_Other,ds,0,0,0);
654: return(0);
655: }
659: /*@C
660: DSCond - Compute the inf-norm condition number of the first matrix
661: as cond(A) = norm(A)*norm(inv(A)).
663: Not Collective
665: Input Parameters:
666: + ds - the direct solver context
667: - cond - the computed condition number
669: Level: advanced
671: .seealso: DSSolve()
672: @*/
673: PetscErrorCode DSCond(DS ds,PetscReal *cond)
674: {
680: if (!ds->ops->cond) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
681: PetscLogEventBegin(DS_Other,ds,0,0,0);
682: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
683: (*ds->ops->cond)(ds,cond);
684: PetscFPTrapPop();
685: PetscLogEventEnd(DS_Other,ds,0,0,0);
686: return(0);
687: }
691: /*@C
692: DSTranslateHarmonic - Computes a translation of the dense system.
694: Logically Collective on DS
696: Input Parameters:
697: + ds - the direct solver context
698: . tau - the translation amount
699: . beta - last component of vector b
700: - recover - boolean flag to indicate whether to recover or not
702: Output Parameters:
703: + g - the computed vector (optional)
704: - gamma - scale factor (optional)
706: Notes:
707: This function is intended for use in the context of Krylov methods only.
708: It computes a translation of a Krylov decomposition in order to extract
709: eigenpair approximations by harmonic Rayleigh-Ritz.
710: The matrix is updated as A + g*b' where g = (A-tau*eye(n))'\b and
711: vector b is assumed to be beta*e_n^T.
713: The gamma factor is defined as sqrt(1+g'*g) and can be interpreted as
714: the factor by which the residual of the Krylov decomposition is scaled.
716: If the recover flag is activated, the computed translation undoes the
717: translation done previously. In that case, parameter tau is ignored.
719: Level: developer
720: @*/
721: PetscErrorCode DSTranslateHarmonic(DS ds,PetscScalar tau,PetscReal beta,PetscBool recover,PetscScalar *g,PetscReal *gamma)
722: {
727: if (!ds->ops->transharm) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
728: PetscLogEventBegin(DS_Other,ds,0,0,0);
729: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
730: (*ds->ops->transharm)(ds,tau,beta,recover,g,gamma);
731: PetscFPTrapPop();
732: PetscLogEventEnd(DS_Other,ds,0,0,0);
733: ds->state = DS_STATE_RAW;
734: PetscObjectStateIncrease((PetscObject)ds);
735: return(0);
736: }
740: /*@C
741: DSTranslateRKS - Computes a modification of the dense system corresponding
742: to an update of the shift in a rational Krylov method.
744: Logically Collective on DS
746: Input Parameters:
747: + ds - the direct solver context
748: - alpha - the translation amount
750: Notes:
751: This function is intended for use in the context of Krylov methods only.
752: It takes the leading (k+1,k) submatrix of A, containing the truncated
753: Rayleigh quotient of a Krylov-Schur relation computed from a shift
754: sigma1 and transforms it to obtain a Krylov relation as if computed
755: from a different shift sigma2. The new matrix is computed as
756: 1.0/alpha*(eye(k)-Q*inv(R)), where [Q,R]=qr(eye(k)-alpha*A) and
757: alpha = sigma1-sigma2.
759: Matrix Q is placed in DS_MAT_Q so that it can be used to update the
760: Krylov basis.
762: Level: developer
763: @*/
764: PetscErrorCode DSTranslateRKS(DS ds,PetscScalar alpha)
765: {
770: if (!ds->ops->transrks) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
771: PetscLogEventBegin(DS_Other,ds,0,0,0);
772: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
773: (*ds->ops->transrks)(ds,alpha);
774: PetscFPTrapPop();
775: PetscLogEventEnd(DS_Other,ds,0,0,0);
776: ds->state = DS_STATE_RAW;
777: ds->compact = PETSC_FALSE;
778: PetscObjectStateIncrease((PetscObject)ds);
779: return(0);
780: }