Actual source code: stsolve.c

  2: /*
  3:     The ST (spectral transformation) interface routines, callable by users.
  4: */

 6:  #include src/st/stimpl.h

 10: /*@
 11:    STApply - Applies the spectral transformation operator to a vector, for
 12:    instance (A - sB)^-1 B in the case of the shift-and-invert tranformation
 13:    and generalized eigenproblem.

 15:    Collective on ST and Vec

 17:    Input Parameters:
 18: +  st - the spectral transformation context
 19: -  x  - input vector

 21:    Output Parameter:
 22: .  y - output vector

 24:    Level: developer

 26: .seealso: STApplyB(), STApplyNoB()
 27: @*/
 28: PetscErrorCode STApply(ST st,Vec x,Vec y)
 29: {

 36:   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");

 38:   if (!st->setupcalled) { STSetUp(st);  }

 40:   PetscLogEventBegin(ST_Apply,st,x,y,0);
 41:   st->applys++;
 42:   (*st->ops->apply)(st,x,y);
 43:   PetscLogEventEnd(ST_Apply,st,x,y,0);
 44:   return(0);
 45: }

 49: /*@
 50:    STApplyB - Applies the B matrix to a vector.

 52:    Collective on ST and Vec

 54:    Input Parameters:
 55: +  st - the spectral transformation context
 56: -  x - input vector

 58:    Output Parameter:
 59: .  y - output vector

 61:    Level: developer

 63: .seealso: STApply(), STApplyNoB()
 64: @*/
 65: PetscErrorCode STApplyB(ST st,Vec x,Vec y)
 66: {

 73:   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");

 75:   if (!st->setupcalled) { STSetUp(st);  }

 77:   if (x->id == st->xid && x->state == st->xstate) {
 78:     VecCopy(st->Bx, y);
 79:     return(0);
 80:   }

 82:   PetscLogEventBegin(ST_ApplyB,st,x,y,0);
 83:   (*st->ops->applyB)(st,x,y);
 84:   PetscLogEventEnd(ST_ApplyB,st,x,y,0);
 85: 
 86:   st->xid = x->id;
 87:   st->xstate = x->state;
 88:   VecCopy(y,st->Bx);
 89:   return(0);
 90: }

 94: /*@
 95:    STApplyNoB - Applies the spectral transformation operator to a vector 
 96:    which has already been multiplied by matrix B. For instance, this routine
 97:    would perform the operation y =(A - sB)^-1 x in the case of the 
 98:    shift-and-invert tranformation and generalized eigenproblem.

100:    Collective on ST and Vec

102:    Input Parameters:
103: +  st - the spectral transformation context
104: -  x  - input vector, where it is assumed that x=Bw for some vector w

106:    Output Parameter:
107: .  y - output vector

109:    Level: developer

111: .seealso: STApply(), STApplyB()
112: @*/
113: PetscErrorCode STApplyNoB(ST st,Vec x,Vec y)
114: {

121:   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");

123:   if (!st->setupcalled) { STSetUp(st);  }

125:   PetscLogEventBegin(ST_ApplyNoB,st,x,y,0);
126:   st->applys++;
127:   (*st->ops->applynoB)(st,x,y);
128:   PetscLogEventEnd(ST_ApplyNoB,st,x,y,0);
129:   return(0);
130: }

134: /*@
135:    STApplyTranspose - Applies the transpose of the operator to a vector, for
136:    instance B^T(A - sB)^-T in the case of the shift-and-invert tranformation
137:    and generalized eigenproblem.

139:    Collective on ST and Vec

141:    Input Parameters:
142: +  st - the spectral transformation context
143: -  x  - input vector

145:    Output Parameter:
146: .  y - output vector

148:    Level: developer

150: .seealso: STApply()
151: @*/
152: PetscErrorCode STApplyTranspose(ST st,Vec x,Vec y)
153: {

160:   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");

162:   if (!st->setupcalled) { STSetUp(st);  }

164:   PetscLogEventBegin(ST_ApplyTranspose,st,x,y,0);
165:   (*st->ops->applytrans)(st,x,y);
166:   PetscLogEventEnd(ST_ApplyTranspose,st,x,y,0);
167:   return(0);
168: }

172: /*@
173:     STComputeExplicitOperator - Computes the explicit operator associated
174:     to the eigenvalue problem with the specified spectral transformation.  

176:     Collective on ST

178:     Input Parameter:
179: .   st - the spectral transform context

181:     Output Parameter:
182: .   mat - the explicit operator

184:     Notes:
185:     This routine builds a matrix containing the explicit operator. For 
186:     example, in generalized problems with shift-and-invert spectral
187:     transformation the result would be matrix (A - s B)^-1 B.
188:     
189:     This computation is done by applying the operator to columns of the 
190:     identity matrix. Note that the result is a dense matrix.

192:     Level: advanced

194: .seealso: STApply()   
195: @*/
196: PetscErrorCode STComputeExplicitOperator(ST st,Mat *mat)
197: {
199:   Vec            in,out;
200:   PetscInt       i,M,m,*rows,start,end;
201:   PetscScalar    *array,one = 1.0;


207:   MatGetVecs(st->A,&in,&out);
208:   VecGetSize(out,&M);
209:   VecGetLocalSize(out,&m);
210:   VecGetOwnershipRange(out,&start,&end);
211:   PetscMalloc(m*sizeof(int),&rows);
212:   for (i=0; i<m; i++) rows[i] = start + i;

214:   MatCreateMPIDense(st->comm,m,m,M,M,PETSC_NULL,mat);

216:   for (i=0; i<M; i++) {
217:     VecSet(in,0.0);
218:     VecSetValues(in,1,&i,&one,INSERT_VALUES);
219:     VecAssemblyBegin(in);
220:     VecAssemblyEnd(in);

222:     STApply(st,in,out);
223: 
224:     VecGetArray(out,&array);
225:     MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
226:     VecRestoreArray(out,&array);
227:   }
228:   PetscFree(rows);
229:   VecDestroy(in);
230:   VecDestroy(out);
231:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
232:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
233:   return(0);
234: }

238: /*@
239:    STNorm - Computes de norm of a vector as the square root of the inner 
240:    product (x,x) as defined by STInnerProduct().

242:    Collective on ST and Vec

244:    Input Parameters:
245: +  st - the spectral transformation context
246: -  x  - input vector

248:    Output Parameter:
249: .  norm - the computed norm

251:    Notes:
252:    This function will usually compute the 2-norm of a vector, ||x||_2. But
253:    this behaviour may be different if using a non-standard inner product changed 
254:    via STSetBilinearForm(). For example, if using the B-inner product for 
255:    positive definite B, (x,y)_B=y^H Bx, then the computed norm is ||x||_B = 
256:    sqrt( x^H Bx ).

258:    Level: developer

260: .seealso: STInnerProduct()
261: @*/
262: PetscErrorCode STNorm(ST st,Vec x,PetscReal *norm)
263: {
265:   PetscScalar    p;

271: 
272:   STInnerProduct(st,x,x,&p);

274:   if (PetscAbsScalar(p)<PETSC_MACHINE_EPSILON)
275:     PetscInfo(st,"Zero norm, either the vector is zero or a semi-inner product is being used\n");

277: #if defined(PETSC_USE_COMPLEX)
278:   if (PetscRealPart(p)<0.0 || PetscAbsReal(PetscImaginaryPart(p))>PETSC_MACHINE_EPSILON)
279:      SETERRQ(1,"STNorm: The inner product is not well defined");
280:   *norm = PetscSqrtScalar(PetscRealPart(p));
281: #else
282:   if (p<0.0) SETERRQ(1,"STNorm: The inner product is not well defined");
283:   *norm = PetscSqrtScalar(p);
284: #endif

286:   return(0);
287: }

291: /*@
292:    STNormBegin - Starts a split phase norm computation.

294:    Input Parameters:
295: +  st   - the spectral transformation context
296: .  x    - input vector
297: -  norm - where the result will go

299:    Level: developer

301:    Notes:
302:    Each call to STNormBegin() should be paired with a call to STNormEnd().

304: .seealso: STNormEnd(), STNorm(), STInnerProduct(), STMInnerProduct(), 
305:           STInnerProductBegin(), STInnerProductEnd()

307: @*/
308: PetscErrorCode STNormBegin(ST st,Vec x,PetscReal *norm)
309: {
311:   PetscScalar    p;

317: 
318:   STInnerProductBegin(st,x,x,&p);

320:   return(0);
321: }

325: /*@
326:    STNormEnd - Ends a split phase norm computation.

328:    Input Parameters:
329: +  st   - the spectral transformation context
330: -  x    - input vector

332:    Output Parameter:
333: .  norm - the computed norm

335:    Level: developer

337:    Notes:
338:    Each call to STNormBegin() should be paired with a call to STNormEnd().

340: .seealso: STNormBegin(), STNorm(), STInnerProduct(), STMInnerProduct(), 
341:           STInnerProductBegin(), STInnerProductEnd()

343: @*/
344: PetscErrorCode STNormEnd(ST st,Vec x,PetscReal *norm)
345: {
347:   PetscScalar    p;

353: 
354:   STInnerProductEnd(st,x,x,&p);

356:   if (PetscAbsScalar(p)<PETSC_MACHINE_EPSILON)
357:     PetscInfo(st,"Zero norm, either the vector is zero or a semi-inner product is being used\n");

359: #if defined(PETSC_USE_COMPLEX)
360:   if (PetscRealPart(p)<0.0 || PetscAbsReal(PetscImaginaryPart(p))>PETSC_MACHINE_EPSILON)
361:      SETERRQ(1,"STNorm: The inner product is not well defined");
362:   *norm = PetscSqrtScalar(PetscRealPart(p));
363: #else
364:   if (p<0.0) SETERRQ(1,"STNorm: The inner product is not well defined");
365:   *norm = PetscSqrtScalar(p);
366: #endif

368:   return(0);
369: }

373: /*@
374:    STInnerProduct - Computes the inner product of two vectors.

376:    Collective on ST and Vec

378:    Input Parameters:
379: +  st - the spectral transformation context
380: .  x  - input vector
381: -  y  - input vector

383:    Output Parameter:
384: .  p - result of the inner product

386:    Notes:
387:    This function will usually compute the standard dot product of vectors
388:    x and y, (x,y)=y^H x. However this behaviour may be different if changed 
389:    via STSetBilinearForm(). This allows use of other inner products such as
390:    the indefinite product y^T x for complex symmetric problems or the
391:    B-inner product for positive definite B, (x,y)_B=y^H Bx.

393:    Level: developer

395: .seealso: STSetBilinearForm(), STApplyB(), VecDot(), STMInnerProduct()
396: @*/
397: PetscErrorCode STInnerProduct(ST st,Vec x,Vec y,PetscScalar *p)
398: {

406: 
407:   PetscLogEventBegin(ST_InnerProduct,st,x,0,0);
408:   st->innerproducts++;
409:   switch (st->bilinear_form) {
410:   case STINNER_HERMITIAN:
411:   case STINNER_SYMMETRIC:
412:     VecCopy(x,st->w);
413:     break;
414:   case STINNER_B_HERMITIAN:
415:   case STINNER_B_SYMMETRIC:
416:     STApplyB(st,x,st->w);
417:     break;
418:   }
419:   switch (st->bilinear_form) {
420:   case STINNER_HERMITIAN:
421:   case STINNER_B_HERMITIAN:
422:     VecDot(st->w,y,p);
423:     break;
424:   case STINNER_SYMMETRIC:
425:   case STINNER_B_SYMMETRIC:
426:     VecTDot(st->w,y,p);
427:     break;
428:   }
429:   PetscLogEventEnd(ST_InnerProduct,st,x,0,0);
430:   return(0);
431: }

435: /*@
436:    STInnerProductBegin - Starts a split phase inner product computation.

438:    Input Parameters:
439: +  st - the spectral transformation context
440: .  x  - the first vector
441: .  y  - the second vector
442: -  p  - where the result will go

444:    Level: developer

446:    Notes:
447:    Each call to STInnerProductBegin() should be paired with a call to STInnerProductEnd().

449: .seealso: STInnerProductEnd(), STInnerProduct(), STNorm(), STNormBegin(), 
450:           STNormEnd(), STMInnerProduct() 

452: @*/
453: PetscErrorCode STInnerProductBegin(ST st,Vec x,Vec y,PetscScalar *p)
454: {

462: 
463:   PetscLogEventBegin(ST_InnerProduct,st,x,0,0);
464:   st->innerproducts++;
465:   switch (st->bilinear_form) {
466:   case STINNER_HERMITIAN:
467:   case STINNER_SYMMETRIC:
468:     VecCopy(x,st->w);
469:     break;
470:   case STINNER_B_HERMITIAN:
471:   case STINNER_B_SYMMETRIC:
472:     STApplyB(st,x,st->w);
473:     break;
474:   }
475:   switch (st->bilinear_form) {
476:   case STINNER_HERMITIAN:
477:   case STINNER_B_HERMITIAN:
478:     VecDotBegin(st->w,y,p);
479:     break;
480:   case STINNER_SYMMETRIC:
481:   case STINNER_B_SYMMETRIC:
482:     VecTDotBegin(st->w,y,p);
483:     break;
484:   }
485:   PetscLogEventEnd(ST_InnerProduct,st,x,0,0);
486:   return(0);
487: }

491: /*@
492:    STInnerProductEnd - Ends a split phase inner product computation.

494:    Input Parameters:
495: +  st - the spectral transformation context
496: .  x  - the first vector
497: -  y  - the second vector

499:    Output Parameter:
500: .  p  - result of the inner product

502:    Level: developer

504:    Notes:
505:    Each call to STInnerProductBegin() should be paired with a call to STInnerProductEnd().

507: .seealso: STInnerProductBegin(), STInnerProduct(), STNorm(), STNormBegin(), 
508:           STNormEnd(), STMInnerProduct() 

510: @*/
511: PetscErrorCode STInnerProductEnd(ST st,Vec x,Vec y,PetscScalar *p)
512: {

520: 
521:   PetscLogEventBegin(ST_InnerProduct,st,x,0,0);
522:   switch (st->bilinear_form) {
523:   case STINNER_HERMITIAN:
524:   case STINNER_B_HERMITIAN:
525:     VecDotEnd(st->w,y,p);
526:     break;
527:   case STINNER_SYMMETRIC:
528:   case STINNER_B_SYMMETRIC:
529:     VecTDotEnd(st->w,y,p);
530:     break;
531:   }
532:   PetscLogEventEnd(ST_InnerProduct,st,x,0,0);
533:   return(0);
534: }

538: /*@
539:    STMInnerProduct - Computes the inner products a vector x with a set of
540:    vectors (columns of Y).

542:    Collective on ST and Vec

544:    Input Parameters:
545: +  st - the spectral transformation context
546: .  n  - number of vectors in y
547: .  x  - the first input vector
548: -  y  - array of vectors

550:    Output Parameter:
551: .  p - result of the inner products

553:    Notes:
554:    This function will usually compute the standard dot product of x and y_i, 
555:    (x,y_i)=y_i^H x, for each column of Y. However this behaviour may be different
556:    if changed via STSetBilinearForm(). This allows use of other inner products 
557:    such as the indefinite product y_i^T x for complex symmetric problems or the
558:    B-inner product for positive definite B, (x,y_i)_B=y_i^H Bx.

560:    Level: developer

562: .seealso: STSetBilinearForm(), STApplyB(), VecMDot(), STInnerProduct()
563: @*/
564: PetscErrorCode STMInnerProduct(ST st,PetscInt n,Vec x,const Vec y[],PetscScalar *p)
565: {

574: 
575:   PetscLogEventBegin(ST_InnerProduct,st,x,0,0);
576:   st->innerproducts += n;
577:   switch (st->bilinear_form) {
578:   case STINNER_HERMITIAN:
579:   case STINNER_SYMMETRIC:
580:     VecCopy(x,st->w);
581:     break;
582:   case STINNER_B_HERMITIAN:
583:   case STINNER_B_SYMMETRIC:
584:     STApplyB(st,x,st->w);
585:     break;
586:   }
587:   switch (st->bilinear_form) {
588:   case STINNER_HERMITIAN:
589:   case STINNER_B_HERMITIAN:
590:     VecMDot(st->w,n,y,p);
591:     break;
592:   case STINNER_SYMMETRIC:
593:   case STINNER_B_SYMMETRIC:
594:     VecMTDot(st->w,n,y,p);
595:     break;
596:   }
597:   PetscLogEventEnd(ST_InnerProduct,st,x,0,0);
598:   return(0);
599: }

603: /*@
604:    STMInnerProductBegin - Starts a split phase multiple inner product computation.

606:    Input Parameters:
607: +  st - the spectral transformation context
608: .  n  - number of vectors in y
609: .  x  - the first input vector
610: .  y  - array of vectors
611: -  p  - where the result will go

613:    Level: developer

615:    Notes:
616:    Each call to STMInnerProductBegin() should be paired with a call to STMInnerProductEnd().

618: .seealso: STMInnerProductEnd(), STMInnerProduct(), STNorm(), STNormBegin(), 
619:           STNormEnd(), STInnerProduct() 

621: @*/
622: PetscErrorCode STMInnerProductBegin(ST st,PetscInt n,Vec x,const Vec y[],PetscScalar *p)
623: {

632: 
633:   PetscLogEventBegin(ST_InnerProduct,st,x,0,0);
634:   st->innerproducts += n;
635:   switch (st->bilinear_form) {
636:   case STINNER_HERMITIAN:
637:   case STINNER_SYMMETRIC:
638:     VecCopy(x,st->w);
639:     break;
640:   case STINNER_B_HERMITIAN:
641:   case STINNER_B_SYMMETRIC:
642:     STApplyB(st,x,st->w);
643:     break;
644:   }
645:   switch (st->bilinear_form) {
646:   case STINNER_HERMITIAN:
647:   case STINNER_B_HERMITIAN:
648:     VecMDotBegin(st->w,n,y,p);
649:     break;
650:   case STINNER_SYMMETRIC:
651:   case STINNER_B_SYMMETRIC:
652:     VecMTDotBegin(st->w,n,y,p);
653:     break;
654:   }
655:   PetscLogEventEnd(ST_InnerProduct,st,x,0,0);
656:   return(0);
657: }

661: /*@
662:    STMInnerProductEnd - Ends a split phase multiple inner product computation.

664:    Input Parameters:
665: +  st - the spectral transformation context
666: .  n  - number of vectors in y
667: .  x  - the first input vector
668: -  y  - array of vectors

670:    Output Parameter:
671: .  p - result of the inner products

673:    Level: developer

675:    Notes:
676:    Each call to STMInnerProductBegin() should be paired with a call to STMInnerProductEnd().

678: .seealso: STMInnerProductBegin(), STMInnerProduct(), STNorm(), STNormBegin(), 
679:           STNormEnd(), STInnerProduct() 

681: @*/
682: PetscErrorCode STMInnerProductEnd(ST st,PetscInt n,Vec x,const Vec y[],PetscScalar *p)
683: {

692: 
693:   PetscLogEventBegin(ST_InnerProduct,st,x,0,0);
694:   switch (st->bilinear_form) {
695:   case STINNER_HERMITIAN:
696:   case STINNER_B_HERMITIAN:
697:     VecMDotEnd(st->w,n,y,p);
698:     break;
699:   case STINNER_SYMMETRIC:
700:   case STINNER_B_SYMMETRIC:
701:     VecMTDotEnd(st->w,n,y,p);
702:     break;
703:   }
704:   PetscLogEventEnd(ST_InnerProduct,st,x,0,0);
705:   return(0);
706: }

710: /*@
711:    STSetUp - Prepares for the use of a spectral transformation.

713:    Collective on ST

715:    Input Parameter:
716: .  st - the spectral transformation context

718:    Level: advanced

720: .seealso: STCreate(), STApply(), STDestroy()
721: @*/
722: PetscErrorCode STSetUp(ST st)
723: {


729:   PetscInfo(st,"Setting up new ST\n");
730:   if (st->setupcalled) return(0);
731:   PetscLogEventBegin(ST_SetUp,st,0,0,0);
732:   if (!st->A) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
733:   if (!st->type_name) {
734:     STSetType(st,STSHIFT);
735:   }
736:   if (st->w) { VecDestroy(st->w); }
737:   if (st->Bx) { VecDestroy(st->Bx); }
738:   MatGetVecs(st->A,&st->w,&st->Bx);
739:   st->xid = 0;
740:   if (st->ops->setup) {
741:     (*st->ops->setup)(st);
742:   }
743:   st->setupcalled = 1;
744:   PetscLogEventEnd(ST_SetUp,st,0,0,0);
745:   return(0);
746: }

750: /*@
751:    STPostSolve - Optional post-solve phase, intended for any actions that must 
752:    be performed on the ST object after the eigensolver has finished.

754:    Collective on ST

756:    Input Parameters:
757: .  st  - the spectral transformation context

759:    Level: developer

761: .seealso: EPSSolve()
762: @*/
763: PetscErrorCode STPostSolve(ST st)
764: {

769:   if (st->ops->postsolve) {
770:     (*st->ops->postsolve)(st);
771:   }

773:   return(0);
774: }

778: /*@
779:    STBackTransform - Back-transformation phase, intended for 
780:    spectral transformations which require to transform the computed 
781:    eigenvalues back to the original eigenvalue problem.

783:    Collective on ST

785:    Input Parameters:
786:    st   - the spectral transformation context
787:    eigr - real part of a computed eigenvalue
788:    eigi - imaginary part of a computed eigenvalue

790:    Level: developer

792: .seealso: EPSBackTransform()
793: @*/
794: PetscErrorCode STBackTransform(ST st,PetscScalar* eigr,PetscScalar* eigi)
795: {

800:   if (st->ops->backtr) {
801:     (*st->ops->backtr)(st,eigr,eigi);
802:   }

804:   return(0);
805: }