Actual source code: ipdot.c
1: /*
2: Dot product routines
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc. See the README file for conditions of use
9: and additional information.
10: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
11: */
13: #include "src/ip/ipimpl.h" /*I "slepcip.h" I*/
17: /*@
18: IPNorm - Computes de norm of a vector as the square root of the inner
19: product (x,x) as defined by IPInnerProduct().
21: Collective on IP and Vec
23: Input Parameters:
24: + ip - the spectral transformation context
25: - x - input vector
27: Output Parameter:
28: . norm - the computed norm
30: Notes:
31: This function will usually compute the 2-norm of a vector, ||x||_2. But
32: this behaviour may be different if using a non-standard inner product changed
33: via IPSetBilinearForm(). For example, if using the B-inner product for
34: positive definite B, (x,y)_B=y^H Bx, then the computed norm is ||x||_B =
35: sqrt( x^H Bx ).
37: Level: developer
39: .seealso: IPInnerProduct()
40: @*/
41: PetscErrorCode IPNorm(IP ip,Vec x,PetscReal *norm)
42: {
44: PetscScalar p;
50:
51: IPInnerProduct(ip,x,x,&p);
53: if (PetscAbsScalar(p)<PETSC_MACHINE_EPSILON)
54: PetscInfo(ip,"Zero norm, either the vector is zero or a semi-inner product is being used\n");
56: #if defined(PETSC_USE_COMPLEX)
57: if (PetscRealPart(p)<0.0 || PetscAbsReal(PetscImaginaryPart(p))>PETSC_MACHINE_EPSILON)
58: SETERRQ(1,"IPNorm: The inner product is not well defined");
59: *norm = PetscSqrtScalar(PetscRealPart(p));
60: #else
61: if (p<0.0) SETERRQ(1,"IPNorm: The inner product is not well defined");
62: *norm = PetscSqrtScalar(p);
63: #endif
65: return(0);
66: }
70: /*@
71: IPNormBegin - Starts a split phase norm computation.
73: Input Parameters:
74: + ip - the spectral transformation context
75: . x - input vector
76: - norm - where the result will go
78: Level: developer
80: Notes:
81: Each call to IPNormBegin() should be paired with a call to IPNormEnd().
83: .seealso: IPNormEnd(), IPNorm(), IPInnerProduct(), IPMInnerProduct(),
84: IPInnerProductBegin(), IPInnerProductEnd()
86: @*/
87: PetscErrorCode IPNormBegin(IP ip,Vec x,PetscReal *norm)
88: {
90: PetscScalar p;
96:
97: IPInnerProductBegin(ip,x,x,&p);
99: return(0);
100: }
104: /*@
105: IPNormEnd - Ends a split phase norm computation.
107: Input Parameters:
108: + ip - the spectral transformation context
109: - x - input vector
111: Output Parameter:
112: . norm - the computed norm
114: Level: developer
116: Notes:
117: Each call to IPNormBegin() should be paired with a call to IPNormEnd().
119: .seealso: IPNormBegin(), IPNorm(), IPInnerProduct(), IPMInnerProduct(),
120: IPInnerProductBegin(), IPInnerProductEnd()
122: @*/
123: PetscErrorCode IPNormEnd(IP ip,Vec x,PetscReal *norm)
124: {
126: PetscScalar p;
132:
133: IPInnerProductEnd(ip,x,x,&p);
135: if (PetscAbsScalar(p)<PETSC_MACHINE_EPSILON)
136: PetscInfo(ip,"Zero norm, either the vector is zero or a semi-inner product is being used\n");
138: #if defined(PETSC_USE_COMPLEX)
139: if (PetscRealPart(p)<0.0 || PetscAbsReal(PetscImaginaryPart(p))>PETSC_MACHINE_EPSILON)
140: SETERRQ(1,"IPNorm: The inner product is not well defined");
141: *norm = PetscSqrtScalar(PetscRealPart(p));
142: #else
143: if (p<0.0) SETERRQ(1,"IPNorm: The inner product is not well defined");
144: *norm = PetscSqrtScalar(p);
145: #endif
147: return(0);
148: }
152: /*@
153: IPInnerProduct - Computes the inner product of two vectors.
155: Collective on IP and Vec
157: Input Parameters:
158: + ip - the spectral transformation context
159: . x - input vector
160: - y - input vector
162: Output Parameter:
163: . p - result of the inner product
165: Notes:
166: This function will usually compute the standard dot product of vectors
167: x and y, (x,y)=y^H x. However this behaviour may be different if changed
168: via IPSetBilinearForm(). This allows use of other inner products such as
169: the indefinite product y^T x for complex symmetric problems or the
170: B-inner product for positive definite B, (x,y)_B=y^H Bx.
172: Level: developer
174: .seealso: IPSetBilinearForm(), IPApplyB(), VecDot(), IPMInnerProduct()
175: @*/
176: PetscErrorCode IPInnerProduct(IP ip,Vec x,Vec y,PetscScalar *p)
177: {
185:
186: PetscLogEventBegin(IP_InnerProduct,ip,x,0,0);
187: ip->innerproducts++;
188: if (ip->matrix) {
189: IPApplyMatrix_Private(ip,x);
190: if (ip->bilinear_form == IPINNER_HERMITIAN) {
191: VecDot(ip->Bx,y,p);
192: } else {
193: VecTDot(ip->Bx,y,p);
194: }
195: } else {
196: if (ip->bilinear_form == IPINNER_HERMITIAN) {
197: VecDot(x,y,p);
198: } else {
199: VecTDot(x,y,p);
200: }
201: }
202: PetscLogEventEnd(IP_InnerProduct,ip,x,0,0);
203: return(0);
204: }
208: /*@
209: IPInnerProductBegin - Starts a split phase inner product computation.
211: Input Parameters:
212: + ip - the spectral transformation context
213: . x - the first vector
214: . y - the second vector
215: - p - where the result will go
217: Level: developer
219: Notes:
220: Each call to IPInnerProductBegin() should be paired with a call to IPInnerProductEnd().
222: .seealso: IPInnerProductEnd(), IPInnerProduct(), IPNorm(), IPNormBegin(),
223: IPNormEnd(), IPMInnerProduct()
225: @*/
226: PetscErrorCode IPInnerProductBegin(IP ip,Vec x,Vec y,PetscScalar *p)
227: {
235:
236: PetscLogEventBegin(IP_InnerProduct,ip,x,0,0);
237: ip->innerproducts++;
238: if (ip->matrix) {
239: IPApplyMatrix_Private(ip,x);
240: if (ip->bilinear_form == IPINNER_HERMITIAN) {
241: VecDotBegin(ip->Bx,y,p);
242: } else {
243: VecTDotBegin(ip->Bx,y,p);
244: }
245: } else {
246: if (ip->bilinear_form == IPINNER_HERMITIAN) {
247: VecDotBegin(x,y,p);
248: } else {
249: VecTDotBegin(x,y,p);
250: }
251: }
252: PetscLogEventEnd(IP_InnerProduct,ip,x,0,0);
253: return(0);
254: }
258: /*@
259: IPInnerProductEnd - Ends a split phase inner product computation.
261: Input Parameters:
262: + ip - the spectral transformation context
263: . x - the first vector
264: - y - the second vector
266: Output Parameter:
267: . p - result of the inner product
269: Level: developer
271: Notes:
272: Each call to IPInnerProductBegin() should be paired with a call to IPInnerProductEnd().
274: .seealso: IPInnerProductBegin(), IPInnerProduct(), IPNorm(), IPNormBegin(),
275: IPNormEnd(), IPMInnerProduct()
277: @*/
278: PetscErrorCode IPInnerProductEnd(IP ip,Vec x,Vec y,PetscScalar *p)
279: {
287:
288: PetscLogEventBegin(IP_InnerProduct,ip,x,0,0);
289: if (ip->matrix) {
290: if (ip->bilinear_form == IPINNER_HERMITIAN) {
291: VecDotEnd(ip->Bx,y,p);
292: } else {
293: VecTDotEnd(ip->Bx,y,p);
294: }
295: } else {
296: if (ip->bilinear_form == IPINNER_HERMITIAN) {
297: VecDotEnd(x,y,p);
298: } else {
299: VecTDotEnd(x,y,p);
300: }
301: }
302: PetscLogEventEnd(IP_InnerProduct,ip,x,0,0);
303: return(0);
304: }
308: /*@
309: IPMInnerProduct - Computes the inner products a vector x with a set of
310: vectors (columns of Y).
312: Collective on IP and Vec
314: Input Parameters:
315: + ip - the spectral transformation context
316: . x - the first input vector
317: . n - number of vectors in y
318: - y - array of vectors
320: Output Parameter:
321: . p - result of the inner products
323: Notes:
324: This function will usually compute the standard dot product of x and y_i,
325: (x,y_i)=y_i^H x, for each column of Y. However this behaviour may be different
326: if changed via IPSetBilinearForm(). This allows use of other inner products
327: such as the indefinite product y_i^T x for complex symmetric problems or the
328: B-inner product for positive definite B, (x,y_i)_B=y_i^H Bx.
330: Level: developer
332: .seealso: IPSetBilinearForm(), IPApplyB(), VecMDot(), IPInnerProduct()
333: @*/
334: PetscErrorCode IPMInnerProduct(IP ip,Vec x,PetscInt n,const Vec y[],PetscScalar *p)
335: {
344:
345: PetscLogEventBegin(IP_InnerProduct,ip,x,0,0);
346: ip->innerproducts += n;
347: if (ip->matrix) {
348: IPApplyMatrix_Private(ip,x);
349: if (ip->bilinear_form == IPINNER_HERMITIAN) {
350: VecMDot(ip->Bx,n,y,p);
351: } else {
352: VecMTDot(ip->Bx,n,y,p);
353: }
354: } else {
355: if (ip->bilinear_form == IPINNER_HERMITIAN) {
356: VecMDot(x,n,y,p);
357: } else {
358: VecMTDot(x,n,y,p);
359: }
360: }
361: PetscLogEventEnd(IP_InnerProduct,ip,x,0,0);
362: return(0);
363: }
367: /*@
368: IPMInnerProductBegin - Starts a split phase multiple inner product computation.
370: Input Parameters:
371: + ip - the spectral transformation context
372: . x - the first input vector
373: . n - number of vectors in y
374: . y - array of vectors
375: - p - where the result will go
377: Level: developer
379: Notes:
380: Each call to IPMInnerProductBegin() should be paired with a call to IPMInnerProductEnd().
382: .seealso: IPMInnerProductEnd(), IPMInnerProduct(), IPNorm(), IPNormBegin(),
383: IPNormEnd(), IPInnerProduct()
385: @*/
386: PetscErrorCode IPMInnerProductBegin(IP ip,Vec x,PetscInt n,const Vec y[],PetscScalar *p)
387: {
396:
397: PetscLogEventBegin(IP_InnerProduct,ip,x,0,0);
398: ip->innerproducts += n;
399: if (ip->matrix) {
400: IPApplyMatrix_Private(ip,x);
401: if (ip->bilinear_form == IPINNER_HERMITIAN) {
402: VecMDotBegin(ip->Bx,n,y,p);
403: } else {
404: VecMTDotBegin(ip->Bx,n,y,p);
405: }
406: } else {
407: if (ip->bilinear_form == IPINNER_HERMITIAN) {
408: VecMDotBegin(x,n,y,p);
409: } else {
410: VecMTDotBegin(x,n,y,p);
411: }
412: }
413: PetscLogEventEnd(IP_InnerProduct,ip,x,0,0);
414: return(0);
415: }
419: /*@
420: IPMInnerProductEnd - Ends a split phase multiple inner product computation.
422: Input Parameters:
423: + ip - the spectral transformation context
424: . x - the first input vector
425: . n - number of vectors in y
426: - y - array of vectors
428: Output Parameter:
429: . p - result of the inner products
431: Level: developer
433: Notes:
434: Each call to IPMInnerProductBegin() should be paired with a call to IPMInnerProductEnd().
436: .seealso: IPMInnerProductBegin(), IPMInnerProduct(), IPNorm(), IPNormBegin(),
437: IPNormEnd(), IPInnerProduct()
439: @*/
440: PetscErrorCode IPMInnerProductEnd(IP ip,Vec x,PetscInt n,const Vec y[],PetscScalar *p)
441: {
450:
451: PetscLogEventBegin(IP_InnerProduct,ip,x,0,0);
452: if (ip->matrix) {
453: if (ip->bilinear_form == IPINNER_HERMITIAN) {
454: VecMDotEnd(ip->Bx,n,y,p);
455: } else {
456: VecMTDotEnd(ip->Bx,n,y,p);
457: }
458: } else {
459: if (ip->bilinear_form == IPINNER_HERMITIAN) {
460: VecMDotEnd(x,n,y,p);
461: } else {
462: VecMTDotEnd(x,n,y,p);
463: }
464: }
465: PetscLogEventEnd(IP_InnerProduct,ip,x,0,0);
466: return(0);
467: }