Actual source code: dsops.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*
  2:    DS operations: DSSolve(), DSVectors(), etc

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, 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 DS

131:    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 DS

232:    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(), DSStateType
243: @*/
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 DS

518:    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(), DSStateType
530: @*/
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: /*@C
554:    DSSort - Sorts the result of DSSolve() according to a given sorting
555:    criterion.

557:    Logically Collective on DS

559:    Input Parameters:
560: +  ds   - the direct solver context
561: .  eigr - array containing the computed eigenvalues (real part)
562: .  eigi - array containing the computed eigenvalues (imaginary part)
563: .  rr   - (optional) array containing auxiliary values (real part)
564: -  ri   - (optional) array containing auxiliary values (imaginary part)

566:    Input/Output Parameter:
567: .  k    - (optional) number of elements in the leading group

569:    Notes:
570:    This routine sorts the arrays provided in eigr and eigi, and also
571:    sorts the dense system stored inside ds (assumed to be in condensed form).
572:    The sorting criterion is specified with DSSetSlepcSC().

574:    If arrays rr and ri are provided, then a (partial) reordering based on these
575:    values rather than on the eigenvalues is performed. In symmetric problems
576:    a total order is obtained (parameter k is ignored), but otherwise the result
577:    is sorted only partially. In this latter case, it is only guaranteed that
578:    all the first k elements satisfy the comparison with any of the last n-k
579:    elements. The output value of parameter k is the final number of elements in
580:    the first set.

582:    Level: intermediate

584: .seealso: DSSolve(), DSSetSlepcSC()
585: @*/
586: PetscErrorCode DSSort(DS ds,PetscScalar *eigr,PetscScalar *eigi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
587: {
589:   PetscInt       i;

593:   DSCheckSolved(ds,1);
596:   if (ds->state==DS_STATE_TRUNCATED) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Cannot sort a truncated DS");
597:   if (!ds->ops->sort) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
598:   if (!ds->sc) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must provide a sorting criterion first");
599:   if (k && !rr) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Argument k can only be used together with rr");

601:   for (i=0;i<ds->n;i++) ds->perm[i] = i;   /* initialize to trivial permutation */
602:   PetscLogEventBegin(DS_Other,ds,0,0,0);
603:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
604:   (*ds->ops->sort)(ds,eigr,eigi,rr,ri,k);
605:   PetscFPTrapPop();
606:   PetscLogEventEnd(DS_Other,ds,0,0,0);
607:   PetscObjectStateIncrease((PetscObject)ds);
608:   return(0);
609: }

613: /*@C
614:    DSVectors - Compute vectors associated to the dense system such
615:    as eigenvectors.

617:    Logically Collective on DS

619:    Input Parameters:
620: +  ds  - the direct solver context
621: -  mat - the matrix, used to indicate which vectors are required

623:    Input/Output Parameter:
624: -  j   - (optional) index of vector to be computed

626:    Output Parameter:
627: .  rnorm - (optional) computed residual norm

629:    Notes:
630:    Allowed values for mat are DS_MAT_X, DS_MAT_Y, DS_MAT_U and DS_MAT_VT, to
631:    compute right or left eigenvectors, or left or right singular vectors,
632:    respectively.

634:    If NULL is passed in argument j then all vectors are computed,
635:    otherwise j indicates which vector must be computed. In real non-symmetric
636:    problems, on exit the index j will be incremented when a complex conjugate
637:    pair is found.

639:    This function can be invoked after the dense problem has been solved,
640:    to get the residual norm estimate of the associated Ritz pair. In that
641:    case, the relevant information is returned in rnorm.

643:    For computing eigenvectors, LAPACK's _trevc is used so the matrix must
644:    be in (quasi-)triangular form, or call DSSolve() first.

646:    Level: intermediate

648: .seealso: DSSolve()
649: @*/
650: PetscErrorCode DSVectors(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
651: {

656:   DSCheckAlloc(ds,1);
658:   if (!ds->ops->vectors) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
659:   if (rnorm && !j) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Must give a value of j");
660:   if (!ds->mat[mat]) { DSAllocateMat_Private(ds,mat); }
661:   PetscLogEventBegin(DS_Vectors,ds,0,0,0);
662:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
663:   (*ds->ops->vectors)(ds,mat,j,rnorm);
664:   PetscFPTrapPop();
665:   PetscLogEventEnd(DS_Vectors,ds,0,0,0);
666:   PetscObjectStateIncrease((PetscObject)ds);
667:   return(0);
668: }

672: /*@
673:    DSNormalize - Normalize a column or all the columns of a matrix. Considers
674:    the case when the columns represent the real and the imaginary part of a vector.

676:    Logically Collective on DS

678:    Input Parameter:
679: +  ds  - the direct solver context
680: .  mat - the matrix to be modified
681: -  col - the column to normalize or -1 to normalize all of them

683:    Notes:
684:    The columns are normalized with respect to the 2-norm.

686:    If col and col+1 (or col-1 and col) represent the real and the imaginary
687:    part of a vector, both columns are scaled.

689:    Level: advanced

691: .seealso: DSMatType
692: @*/
693: PetscErrorCode DSNormalize(DS ds,DSMatType mat,PetscInt col)
694: {

699:   DSCheckSolved(ds,1);
702:   if (!ds->ops->normalize) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
703:   if (col<-1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"col should be at least minus one");
704:   PetscLogEventBegin(DS_Other,ds,0,0,0);
705:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
706:   (*ds->ops->normalize)(ds,mat,col);
707:   PetscFPTrapPop();
708:   PetscLogEventEnd(DS_Other,ds,0,0,0);
709:   return(0);
710: }

714: /*@
715:    DSUpdateExtraRow - Performs all necessary operations so that the extra
716:    row gets up-to-date after a call to DSSolve().

718:    Not Collective

720:    Input Parameters:
721: .  ds - the direct solver context

723:    Level: advanced

725: .seealso: DSSolve(), DSSetExtraRow()
726: @*/
727: PetscErrorCode DSUpdateExtraRow(DS ds)
728: {

733:   DSCheckAlloc(ds,1);
734:   if (!ds->ops->update) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
735:   if (!ds->extrarow) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Should have called DSSetExtraRow");
736:   PetscLogEventBegin(DS_Other,ds,0,0,0);
737:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
738:   (*ds->ops->update)(ds);
739:   PetscFPTrapPop();
740:   PetscLogEventEnd(DS_Other,ds,0,0,0);
741:   return(0);
742: }

746: /*@
747:    DSCond - Compute the inf-norm condition number of the first matrix
748:    as cond(A) = norm(A)*norm(inv(A)).

750:    Not Collective

752:    Input Parameters:
753: +  ds - the direct solver context
754: -  cond - the computed condition number

756:    Level: advanced

758: .seealso: DSSolve()
759: @*/
760: PetscErrorCode DSCond(DS ds,PetscReal *cond)
761: {

766:   DSCheckAlloc(ds,1);
768:   if (!ds->ops->cond) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
769:   PetscLogEventBegin(DS_Other,ds,0,0,0);
770:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
771:   (*ds->ops->cond)(ds,cond);
772:   PetscFPTrapPop();
773:   PetscLogEventEnd(DS_Other,ds,0,0,0);
774:   return(0);
775: }

779: /*@C
780:    DSTranslateHarmonic - Computes a translation of the dense system.

782:    Logically Collective on DS

784:    Input Parameters:
785: +  ds      - the direct solver context
786: .  tau     - the translation amount
787: .  beta    - last component of vector b
788: -  recover - boolean flag to indicate whether to recover or not

790:    Output Parameters:
791: +  g       - the computed vector (optional)
792: -  gamma   - scale factor (optional)

794:    Notes:
795:    This function is intended for use in the context of Krylov methods only.
796:    It computes a translation of a Krylov decomposition in order to extract
797:    eigenpair approximations by harmonic Rayleigh-Ritz.
798:    The matrix is updated as A + g*b' where g = (A-tau*eye(n))'\b and
799:    vector b is assumed to be beta*e_n^T.

801:    The gamma factor is defined as sqrt(1+g'*g) and can be interpreted as
802:    the factor by which the residual of the Krylov decomposition is scaled.

804:    If the recover flag is activated, the computed translation undoes the
805:    translation done previously. In that case, parameter tau is ignored.

807:    Level: developer
808: @*/
809: PetscErrorCode DSTranslateHarmonic(DS ds,PetscScalar tau,PetscReal beta,PetscBool recover,PetscScalar *g,PetscReal *gamma)
810: {

815:   DSCheckAlloc(ds,1);
816:   if (!ds->ops->transharm) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
817:   PetscLogEventBegin(DS_Other,ds,0,0,0);
818:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
819:   (*ds->ops->transharm)(ds,tau,beta,recover,g,gamma);
820:   PetscFPTrapPop();
821:   PetscLogEventEnd(DS_Other,ds,0,0,0);
822:   ds->state = DS_STATE_RAW;
823:   PetscObjectStateIncrease((PetscObject)ds);
824:   return(0);
825: }

829: /*@
830:    DSTranslateRKS - Computes a modification of the dense system corresponding
831:    to an update of the shift in a rational Krylov method.

833:    Logically Collective on DS

835:    Input Parameters:
836: +  ds    - the direct solver context
837: -  alpha - the translation amount

839:    Notes:
840:    This function is intended for use in the context of Krylov methods only.
841:    It takes the leading (k+1,k) submatrix of A, containing the truncated
842:    Rayleigh quotient of a Krylov-Schur relation computed from a shift
843:    sigma1 and transforms it to obtain a Krylov relation as if computed
844:    from a different shift sigma2. The new matrix is computed as
845:    1.0/alpha*(eye(k)-Q*inv(R)), where [Q,R]=qr(eye(k)-alpha*A) and
846:    alpha = sigma1-sigma2.

848:    Matrix Q is placed in DS_MAT_Q so that it can be used to update the
849:    Krylov basis.

851:    Level: developer
852: @*/
853: PetscErrorCode DSTranslateRKS(DS ds,PetscScalar alpha)
854: {

859:   DSCheckAlloc(ds,1);
860:   if (!ds->ops->transrks) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DS type %s",((PetscObject)ds)->type_name);
861:   PetscLogEventBegin(DS_Other,ds,0,0,0);
862:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
863:   (*ds->ops->transrks)(ds,alpha);
864:   PetscFPTrapPop();
865:   PetscLogEventEnd(DS_Other,ds,0,0,0);
866:   ds->state   = DS_STATE_RAW;
867:   ds->compact = PETSC_FALSE;
868:   PetscObjectStateIncrease((PetscObject)ds);
869:   return(0);
870: }

874: /*@
875:    DSCopyMat - Copies the contents of a sequential dense Mat object to
876:    the indicated DS matrix, or vice versa. 

878:    Not Collective

880:    Input Parameters:
881: +  ds   - the direct solver context
882: .  m    - the requested matrix
883: .  mr   - first row of m to be considered
884: .  mc   - first column of m to be considered
885: .  A    - Mat object
886: .  Ar   - first row of A to be considered
887: .  Ac   - first column of A to be considered
888: .  rows - number of rows to copy
889: .  cols - number of columns to copy
890: -  out  - whether the data is copied out of the DS

892:    Note:
893:    If out=true, the values of the DS matrix m are copied to A, otherwise
894:    the entries of A are copied to the DS.

896:    Level: developer

898: .seealso: DSGetMat()
899: @*/
900: PetscErrorCode DSCopyMat(DS ds,DSMatType m,PetscInt mr,PetscInt mc,Mat A,PetscInt Ar,PetscInt Ac,PetscInt rows,PetscInt cols,PetscBool out)
901: {
903:   PetscInt       j,mrows,mcols,arows,acols;
904:   PetscScalar    *pA,*M;

908:   DSCheckAlloc(ds,1);
918:   if (m>=DS_NUM_MAT) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid matrix");
919:   if (!ds->mat[m]) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Requested matrix was not created in this DS");
920:   if (!rows || !cols) return(0);

922:   mrows = PetscMax(ds->n,ds->t);
923:   mcols = ds->m? ds->m: ds->n;
924:   MatGetSize(A,&arows,&acols);
925:   if (mr<0 || mr>=mrows) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial row in m");
926:   if (mc<0 || mc>=mcols) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial column in m");
927:   if (Ar<0 || Ar>=arows) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial row in A");
928:   if (Ac<0 || Ac>=acols) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial column in A");
929:   if (mr+rows>mrows || Ar+rows>arows) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid number of rows");
930:   if (mc+cols>mcols || Ac+cols>acols) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid number of columns");

932:   M  = ds->mat[m];
933:   MatDenseGetArray(A,&pA);
934:   for (j=0;j<cols;j++) {
935:     if (out) {
936:       PetscMemcpy(pA+(Ac+j)*arows+Ar,M+(mc+j)*ds->ld+mr,rows*sizeof(PetscScalar));
937:     } else {
938:       PetscMemcpy(M+(mc+j)*ds->ld+mr,pA+(Ac+j)*arows+Ar,rows*sizeof(PetscScalar));
939:     }
940:   }
941:   MatDenseRestoreArray(A,&pA);
942:   return(0);
943: }