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->ops->apply)(st,x,y);
42: PetscLogEventEnd(ST_Apply,st,x,y,0);
43: return(0);
44: }
48: /*@
49: STApplyB - Applies the B matrix to a vector.
51: Collective on ST and Vec
53: Input Parameters:
54: + st - the spectral transformation context
55: - x - input vector
57: Output Parameter:
58: . y - output vector
60: Level: developer
62: .seealso: STApply(), STApplyNoB()
63: @*/
64: PetscErrorCode STApplyB(ST st,Vec x,Vec y)
65: {
72: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
74: if (!st->setupcalled) { STSetUp(st); }
76: if (x->id == st->xid && x->state == st->xstate) {
77: VecCopy(st->Bx, y);
78: return(0);
79: }
81: PetscLogEventBegin(ST_ApplyB,st,x,y,0);
82: (*st->ops->applyB)(st,x,y);
83: PetscLogEventEnd(ST_ApplyB,st,x,y,0);
84:
85: st->xid = x->id;
86: st->xstate = x->state;
87: VecCopy(y,st->Bx);
88: return(0);
89: }
93: /*@
94: STApplyNoB - Applies the spectral transformation operator to a vector
95: which has already been multiplied by matrix B. For instance, this routine
96: would perform the operation y =(A - sB)^-1 x in the case of the
97: shift-and-invert tranformation and generalized eigenproblem.
99: Collective on ST and Vec
101: Input Parameters:
102: + st - the spectral transformation context
103: - x - input vector, where it is assumed that x=Bw for some vector w
105: Output Parameter:
106: . y - output vector
108: Level: developer
110: .seealso: STApply(), STApplyB()
111: @*/
112: PetscErrorCode STApplyNoB(ST st,Vec x,Vec y)
113: {
120: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
122: if (!st->setupcalled) { STSetUp(st); }
124: PetscLogEventBegin(ST_ApplyNoB,st,x,y,0);
125: (*st->ops->applynoB)(st,x,y);
126: PetscLogEventEnd(ST_ApplyNoB,st,x,y,0);
127: return(0);
128: }
132: /*@
133: STApplyTranspose - Applies the transpose of the operator to a vector, for
134: instance B^T(A - sB)^-T in the case of the shift-and-invert tranformation
135: and generalized eigenproblem.
137: Collective on ST and Vec
139: Input Parameters:
140: + st - the spectral transformation context
141: - x - input vector
143: Output Parameter:
144: . y - output vector
146: Level: developer
148: .seealso: STApply()
149: @*/
150: PetscErrorCode STApplyTranspose(ST st,Vec x,Vec y)
151: {
158: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
160: if (!st->setupcalled) { STSetUp(st); }
162: PetscLogEventBegin(ST_ApplyTranspose,st,x,y,0);
163: (*st->ops->applytrans)(st,x,y);
164: PetscLogEventEnd(ST_ApplyTranspose,st,x,y,0);
165: return(0);
166: }
170: /*@
171: STComputeExplicitOperator - Computes the explicit operator associated
172: to the eigenvalue problem with the specified spectral transformation.
174: Collective on ST
176: Input Parameter:
177: . st - the spectral transform context
179: Output Parameter:
180: . mat - the explicit operator
182: Notes:
183: This routine builds a matrix containing the explicit operator. For
184: example, in generalized problems with shift-and-invert spectral
185: transformation the result would be matrix (A - s B)^-1 B.
186:
187: This computation is done by applying the operator to columns of the
188: identity matrix. Note that the result is a dense matrix.
190: Level: advanced
192: .seealso: STApply()
193: @*/
194: PetscErrorCode STComputeExplicitOperator(ST st,Mat *mat)
195: {
197: Vec in,out;
198: PetscInt i,M,m,*rows,start,end;
199: PetscScalar *array,one = 1.0;
205: MatGetVecs(st->A,&in,&out);
206: VecGetSize(out,&M);
207: VecGetLocalSize(out,&m);
208: VecGetOwnershipRange(out,&start,&end);
209: PetscMalloc(m*sizeof(int),&rows);
210: for (i=0; i<m; i++) rows[i] = start + i;
212: MatCreateMPIDense(st->comm,m,m,M,M,PETSC_NULL,mat);
214: for (i=0; i<M; i++) {
215: VecSet(in,0.0);
216: VecSetValues(in,1,&i,&one,INSERT_VALUES);
217: VecAssemblyBegin(in);
218: VecAssemblyEnd(in);
220: STApply(st,in,out);
221:
222: VecGetArray(out,&array);
223: MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
224: VecRestoreArray(out,&array);
225: }
226: PetscFree(rows);
227: VecDestroy(in);
228: VecDestroy(out);
229: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
230: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
231: return(0);
232: }
236: /*@
237: STNorm - Computes de norm of a vector as the square root of the inner
238: product (x,x) as defined by STInnerProduct().
240: Collective on ST and Vec
242: Input Parameters:
243: + st - the spectral transformation context
244: - x - input vector
246: Output Parameter:
247: . norm - the computed norm
249: Notes:
250: This function will usually compute the 2-norm of a vector, ||x||_2. But
251: this behaviour may be different if using a non-standard inner product changed
252: via STSetBilinearForm(). For example, if using the B-inner product for
253: positive definite B, (x,y)_B=y^H Bx, then the computed norm is ||x||_B =
254: sqrt( x^H Bx ).
256: Level: developer
258: .seealso: STInnerProduct()
259: @*/
260: PetscErrorCode STNorm(ST st,Vec x,PetscReal *norm)
261: {
263: PetscScalar p;
269:
270: STInnerProduct(st,x,x,&p);
272: if (PetscAbsScalar(p)<PETSC_MACHINE_EPSILON)
273: PetscInfo(st,"Zero norm, either the vector is zero or a semi-inner product is being used\n");
275: #if defined(PETSC_USE_COMPLEX)
276: if (PetscRealPart(p)<0.0 || PetscAbsReal(PetscImaginaryPart(p))>PETSC_MACHINE_EPSILON)
277: SETERRQ(1,"STNorm: The inner product is not well defined");
278: *norm = PetscSqrtScalar(PetscRealPart(p));
279: #else
280: if (p<0.0) SETERRQ(1,"STNorm: The inner product is not well defined");
281: *norm = PetscSqrtScalar(p);
282: #endif
284: return(0);
285: }
289: /*@
290: STNormBegin - Starts a split phase norm computation.
292: Input Parameters:
293: + st - the spectral transformation context
294: . x - input vector
295: - norm - where the result will go
297: Level: developer
299: Notes:
300: Each call to STNormBegin() should be paired with a call to STNormEnd().
302: .seealso: STNormEnd(), STNorm(), STInnerProduct(), STMInnerProduct(),
303: STInnerProductBegin(), STInnerProductEnd()
305: @*/
306: PetscErrorCode STNormBegin(ST st,Vec x,PetscReal *norm)
307: {
309: PetscScalar p;
315:
316: STInnerProductBegin(st,x,x,&p);
318: return(0);
319: }
323: /*@
324: STNormEnd - Ends a split phase norm computation.
326: Input Parameters:
327: + st - the spectral transformation context
328: - x - input vector
330: Output Parameter:
331: . norm - the computed norm
333: Level: developer
335: Notes:
336: Each call to STNormBegin() should be paired with a call to STNormEnd().
338: .seealso: STNormBegin(), STNorm(), STInnerProduct(), STMInnerProduct(),
339: STInnerProductBegin(), STInnerProductEnd()
341: @*/
342: PetscErrorCode STNormEnd(ST st,Vec x,PetscReal *norm)
343: {
345: PetscScalar p;
351:
352: STInnerProductEnd(st,x,x,&p);
354: if (PetscAbsScalar(p)<PETSC_MACHINE_EPSILON)
355: PetscInfo(st,"Zero norm, either the vector is zero or a semi-inner product is being used\n");
357: #if defined(PETSC_USE_COMPLEX)
358: if (PetscRealPart(p)<0.0 || PetscAbsReal(PetscImaginaryPart(p))>PETSC_MACHINE_EPSILON)
359: SETERRQ(1,"STNorm: The inner product is not well defined");
360: *norm = PetscSqrtScalar(PetscRealPart(p));
361: #else
362: if (p<0.0) SETERRQ(1,"STNorm: The inner product is not well defined");
363: *norm = PetscSqrtScalar(p);
364: #endif
366: return(0);
367: }
371: /*@
372: STInnerProduct - Computes the inner product of two vectors.
374: Collective on ST and Vec
376: Input Parameters:
377: + st - the spectral transformation context
378: . x - input vector
379: - y - input vector
381: Output Parameter:
382: . p - result of the inner product
384: Notes:
385: This function will usually compute the standard dot product of vectors
386: x and y, (x,y)=y^H x. However this behaviour may be different if changed
387: via STSetBilinearForm(). This allows use of other inner products such as
388: the indefinite product y^T x for complex symmetric problems or the
389: B-inner product for positive definite B, (x,y)_B=y^H Bx.
391: Level: developer
393: .seealso: STSetBilinearForm(), STApplyB(), VecDot(), STMInnerProduct()
394: @*/
395: PetscErrorCode STInnerProduct(ST st,Vec x,Vec y,PetscScalar *p)
396: {
404:
405: PetscLogEventBegin(ST_InnerProduct,st,x,0,0);
406: switch (st->bilinear_form) {
407: case STINNER_HERMITIAN:
408: case STINNER_SYMMETRIC:
409: VecCopy(x,st->w);
410: break;
411: case STINNER_B_HERMITIAN:
412: case STINNER_B_SYMMETRIC:
413: STApplyB(st,x,st->w);
414: break;
415: }
416: switch (st->bilinear_form) {
417: case STINNER_HERMITIAN:
418: case STINNER_B_HERMITIAN:
419: VecDot(st->w,y,p);
420: break;
421: case STINNER_SYMMETRIC:
422: case STINNER_B_SYMMETRIC:
423: VecTDot(st->w,y,p);
424: break;
425: }
426: PetscLogEventEnd(ST_InnerProduct,st,x,0,0);
427: return(0);
428: }
432: /*@
433: STInnerProductBegin - Starts a split phase inner product computation.
435: Input Parameters:
436: + st - the spectral transformation context
437: . x - the first vector
438: . y - the second vector
439: - p - where the result will go
441: Level: developer
443: Notes:
444: Each call to STInnerProductBegin() should be paired with a call to STInnerProductEnd().
446: .seealso: STInnerProductEnd(), STInnerProduct(), STNorm(), STNormBegin(),
447: STNormEnd(), STMInnerProduct()
449: @*/
450: PetscErrorCode STInnerProductBegin(ST st,Vec x,Vec y,PetscScalar *p)
451: {
459:
460: PetscLogEventBegin(ST_InnerProduct,st,x,0,0);
461: switch (st->bilinear_form) {
462: case STINNER_HERMITIAN:
463: case STINNER_SYMMETRIC:
464: VecCopy(x,st->w);
465: break;
466: case STINNER_B_HERMITIAN:
467: case STINNER_B_SYMMETRIC:
468: STApplyB(st,x,st->w);
469: break;
470: }
471: switch (st->bilinear_form) {
472: case STINNER_HERMITIAN:
473: case STINNER_B_HERMITIAN:
474: VecDotBegin(st->w,y,p);
475: break;
476: case STINNER_SYMMETRIC:
477: case STINNER_B_SYMMETRIC:
478: VecTDotBegin(st->w,y,p);
479: break;
480: }
481: PetscLogEventEnd(ST_InnerProduct,st,x,0,0);
482: return(0);
483: }
487: /*@
488: STInnerProductEnd - Ends a split phase inner product computation.
490: Input Parameters:
491: + st - the spectral transformation context
492: . x - the first vector
493: - y - the second vector
495: Output Parameter:
496: . p - result of the inner product
498: Level: developer
500: Notes:
501: Each call to STInnerProductBegin() should be paired with a call to STInnerProductEnd().
503: .seealso: STInnerProductBegin(), STInnerProduct(), STNorm(), STNormBegin(),
504: STNormEnd(), STMInnerProduct()
506: @*/
507: PetscErrorCode STInnerProductEnd(ST st,Vec x,Vec y,PetscScalar *p)
508: {
516:
517: PetscLogEventBegin(ST_InnerProduct,st,x,0,0);
518: switch (st->bilinear_form) {
519: case STINNER_HERMITIAN:
520: case STINNER_B_HERMITIAN:
521: VecDotEnd(st->w,y,p);
522: break;
523: case STINNER_SYMMETRIC:
524: case STINNER_B_SYMMETRIC:
525: VecTDotEnd(st->w,y,p);
526: break;
527: }
528: PetscLogEventEnd(ST_InnerProduct,st,x,0,0);
529: return(0);
530: }
534: /*@
535: STMInnerProduct - Computes the inner products a vector x with a set of
536: vectors (columns of Y).
538: Collective on ST and Vec
540: Input Parameters:
541: + st - the spectral transformation context
542: . n - number of vectors in y
543: . x - the first input vector
544: - y - array of vectors
546: Output Parameter:
547: . p - result of the inner products
549: Notes:
550: This function will usually compute the standard dot product of x and y_i,
551: (x,y_i)=y_i^H x, for each column of Y. However this behaviour may be different
552: if changed via STSetBilinearForm(). This allows use of other inner products
553: such as the indefinite product y_i^T x for complex symmetric problems or the
554: B-inner product for positive definite B, (x,y_i)_B=y_i^H Bx.
556: Level: developer
558: .seealso: STSetBilinearForm(), STApplyB(), VecMDot(), STInnerProduct()
559: @*/
560: PetscErrorCode STMInnerProduct(ST st,PetscInt n,Vec x,const Vec y[],PetscScalar *p)
561: {
570:
571: PetscLogEventBegin(ST_InnerProduct,st,x,0,0);
572: switch (st->bilinear_form) {
573: case STINNER_HERMITIAN:
574: case STINNER_SYMMETRIC:
575: VecCopy(x,st->w);
576: break;
577: case STINNER_B_HERMITIAN:
578: case STINNER_B_SYMMETRIC:
579: STApplyB(st,x,st->w);
580: break;
581: }
582: switch (st->bilinear_form) {
583: case STINNER_HERMITIAN:
584: case STINNER_B_HERMITIAN:
585: VecMDot(st->w,n,y,p);
586: break;
587: case STINNER_SYMMETRIC:
588: case STINNER_B_SYMMETRIC:
589: VecMTDot(st->w,n,y,p);
590: break;
591: }
592: PetscLogEventEnd(ST_InnerProduct,st,x,0,0);
593: return(0);
594: }
598: /*@
599: STMInnerProductBegin - Starts a split phase multiple inner product computation.
601: Input Parameters:
602: + st - the spectral transformation context
603: . n - number of vectors in y
604: . x - the first input vector
605: . y - array of vectors
606: - p - where the result will go
608: Level: developer
610: Notes:
611: Each call to STMInnerProductBegin() should be paired with a call to STMInnerProductEnd().
613: .seealso: STMInnerProductEnd(), STMInnerProduct(), STNorm(), STNormBegin(),
614: STNormEnd(), STInnerProduct()
616: @*/
617: PetscErrorCode STMInnerProductBegin(ST st,PetscInt n,Vec x,const Vec y[],PetscScalar *p)
618: {
620: int i;
621: PetscTruth mdot;
629:
630: PetscLogEventBegin(ST_InnerProduct,st,x,0,0);
631: switch (st->bilinear_form) {
632: case STINNER_HERMITIAN:
633: case STINNER_SYMMETRIC:
634: VecCopy(x,st->w);
635: break;
636: case STINNER_B_HERMITIAN:
637: case STINNER_B_SYMMETRIC:
638: STApplyB(st,x,st->w);
639: break;
640: }
641: switch (st->bilinear_form) {
642: case STINNER_HERMITIAN:
643: case STINNER_B_HERMITIAN:
644: PetscOptionsHasName(st->prefix,"-mdot",&mdot);
645: if (mdot) {
646: VecMDotBegin(st->w,n,y,p);
647: } else {
648: for (i=0;i<n;i++) {
649: VecDotBegin(st->w,y[i],p+i);
650: }
651: }
652: break;
653: case STINNER_SYMMETRIC:
654: case STINNER_B_SYMMETRIC:
655: PetscOptionsHasName(st->prefix,"-mdot",&mdot);
656: if (mdot) {
657: VecMTDotBegin(st->w,n,y,p);
658: } else {
659: for (i=0;i<n;i++) {
660: VecTDotBegin(st->w,y[i],p+i);
661: }
662: }
663: break;
664: }
665: PetscLogEventEnd(ST_InnerProduct,st,x,0,0);
666: return(0);
667: }
671: /*@
672: STMInnerProductEnd - Ends a split phase multiple inner product computation.
674: Input Parameters:
675: + st - the spectral transformation context
676: . n - number of vectors in y
677: . x - the first input vector
678: - y - array of vectors
680: Output Parameter:
681: . p - result of the inner products
683: Level: developer
685: Notes:
686: Each call to STMInnerProductBegin() should be paired with a call to STMInnerProductEnd().
688: .seealso: STMInnerProductBegin(), STMInnerProduct(), STNorm(), STNormBegin(),
689: STNormEnd(), STInnerProduct()
691: @*/
692: PetscErrorCode STMInnerProductEnd(ST st,PetscInt n,Vec x,const Vec y[],PetscScalar *p)
693: {
695: int i;
696: PetscTruth mdot;
704:
705: PetscLogEventBegin(ST_InnerProduct,st,x,0,0);
706: switch (st->bilinear_form) {
707: case STINNER_HERMITIAN:
708: case STINNER_B_HERMITIAN:
709: PetscOptionsHasName(st->prefix,"-mdot",&mdot);
710: if (mdot) {
711: VecMDotEnd(st->w,n,y,p);
712: } else {
713: for (i=0;i<n;i++) {
714: VecDotEnd(st->w,y[i],p+i);
715: }
716: }
717: break;
718: case STINNER_SYMMETRIC:
719: case STINNER_B_SYMMETRIC:
720: PetscOptionsHasName(st->prefix,"-mdot",&mdot);
721: if (mdot) {
722: VecMTDotEnd(st->w,n,y,p);
723: } else {
724: for (i=0;i<n;i++) {
725: VecTDotEnd(st->w,y[i],p+i);
726: }
727: }
728: break;
729: }
730: PetscLogEventEnd(ST_InnerProduct,st,x,0,0);
731: return(0);
732: }
736: /*@
737: STSetUp - Prepares for the use of a spectral transformation.
739: Collective on ST
741: Input Parameter:
742: . st - the spectral transformation context
744: Level: advanced
746: .seealso: STCreate(), STApply(), STDestroy()
747: @*/
748: PetscErrorCode STSetUp(ST st)
749: {
755: PetscInfo(st,"Setting up new ST\n");
756: if (st->setupcalled) return(0);
757: PetscLogEventBegin(ST_SetUp,st,0,0,0);
758: if (!st->A) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
759: if (!st->type_name) {
760: STSetType(st,STSHIFT);
761: }
762: if (st->w) { VecDestroy(st->w); }
763: if (st->Bx) { VecDestroy(st->Bx); }
764: MatGetVecs(st->A,&st->w,&st->Bx);
765: st->xid = 0;
766: if (st->ops->setup) {
767: (*st->ops->setup)(st);
768: }
769: st->setupcalled = 1;
770: PetscLogEventEnd(ST_SetUp,st,0,0,0);
771: return(0);
772: }
776: /*@
777: STPostSolve - Optional post-solve phase, intended for any actions that must
778: be performed on the ST object after the eigensolver has finished.
780: Collective on ST
782: Input Parameters:
783: . st - the spectral transformation context
785: Level: developer
787: .seealso: EPSSolve()
788: @*/
789: PetscErrorCode STPostSolve(ST st)
790: {
795: if (st->ops->postsolve) {
796: (*st->ops->postsolve)(st);
797: }
799: return(0);
800: }
804: /*@
805: STBackTransform - Back-transformation phase, intended for
806: spectral transformations which require to transform the computed
807: eigenvalues back to the original eigenvalue problem.
809: Collective on ST
811: Input Parameters:
812: st - the spectral transformation context
813: eigr - real part of a computed eigenvalue
814: eigi - imaginary part of a computed eigenvalue
816: Level: developer
818: .seealso: EPSBackTransform()
819: @*/
820: PetscErrorCode STBackTransform(ST st,PetscScalar* eigr,PetscScalar* eigi)
821: {
826: if (st->ops->backtr) {
827: (*st->ops->backtr)(st,eigr,eigi);
828: }
830: return(0);
831: }