Actual source code: ipdot.c
1: /*
2: Dot product routines
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
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 private/ipimpl.h
28: /*@
29: IPNorm - Computes de norm of a vector as the square root of the inner
30: product (x,x) as defined by IPInnerProduct().
32: Collective on IP and Vec
34: Input Parameters:
35: + ip - the inner product context
36: - x - input vector
38: Output Parameter:
39: . norm - the computed norm
41: Notes:
42: This function will usually compute the 2-norm of a vector, ||x||_2. But
43: this behaviour may be different if using a non-standard inner product changed
44: via IPSetBilinearForm(). For example, if using the B-inner product for
45: positive definite B, (x,y)_B=y^H Bx, then the computed norm is ||x||_B =
46: sqrt( x^H Bx ).
48: Level: developer
50: .seealso: IPInnerProduct()
51: @*/
52: PetscErrorCode IPNorm(IP ip,Vec x,PetscReal *norm)
53: {
55: PetscScalar p;
61:
62: if (!ip->matrix && ip->bilinear_form == IPINNER_HERMITIAN) {
63: VecNorm(x,NORM_2,norm);
64: } else {
65: IPInnerProduct(ip,x,x,&p);
66: if (PetscAbsScalar(p)<PETSC_MACHINE_EPSILON)
67: PetscInfo(ip,"Zero norm, either the vector is zero or a semi-inner product is being used\n");
68: #if defined(PETSC_USE_COMPLEX)
69: if (PetscRealPart(p)<0.0 || PetscAbsReal(PetscImaginaryPart(p))>PETSC_MACHINE_EPSILON)
70: SETERRQ(1,"IPNorm: The inner product is not well defined");
71: *norm = PetscSqrtScalar(PetscRealPart(p));
72: #else
73: if (p<0.0) SETERRQ(1,"IPNorm: The inner product is not well defined");
74: *norm = PetscSqrtScalar(p);
75: #endif
76: }
78: return(0);
79: }
83: /*@
84: IPNormBegin - Starts a split phase norm computation.
86: Input Parameters:
87: + ip - the inner product context
88: . x - input vector
89: - norm - where the result will go
91: Level: developer
93: Notes:
94: Each call to IPNormBegin() should be paired with a call to IPNormEnd().
96: .seealso: IPNormEnd(), IPNorm(), IPInnerProduct(), IPMInnerProduct(),
97: IPInnerProductBegin(), IPInnerProductEnd()
99: @*/
100: PetscErrorCode IPNormBegin(IP ip,Vec x,PetscReal *norm)
101: {
103: PetscScalar p;
109:
110: if (!ip->matrix && ip->bilinear_form == IPINNER_HERMITIAN) {
111: VecNormBegin(x,NORM_2,norm);
112: } else {
113: IPInnerProductBegin(ip,x,x,&p);
114: }
116: return(0);
117: }
121: /*@
122: IPNormEnd - Ends a split phase norm computation.
124: Input Parameters:
125: + ip - the inner product context
126: - x - input vector
128: Output Parameter:
129: . norm - the computed norm
131: Level: developer
133: Notes:
134: Each call to IPNormBegin() should be paired with a call to IPNormEnd().
136: .seealso: IPNormBegin(), IPNorm(), IPInnerProduct(), IPMInnerProduct(),
137: IPInnerProductBegin(), IPInnerProductEnd()
139: @*/
140: PetscErrorCode IPNormEnd(IP ip,Vec x,PetscReal *norm)
141: {
143: PetscScalar p;
149:
150: if (!ip->matrix && ip->bilinear_form == IPINNER_HERMITIAN) {
151: VecNormEnd(x,NORM_2,norm);
152: } else {
153: IPInnerProductEnd(ip,x,x,&p);
154: if (PetscAbsScalar(p)<PETSC_MACHINE_EPSILON)
155: PetscInfo(ip,"Zero norm, either the vector is zero or a semi-inner product is being used\n");
157: #if defined(PETSC_USE_COMPLEX)
158: if (PetscRealPart(p)<0.0 || PetscAbsReal(PetscImaginaryPart(p))>PETSC_MACHINE_EPSILON)
159: SETERRQ(1,"IPNorm: The inner product is not well defined");
160: *norm = PetscSqrtScalar(PetscRealPart(p));
161: #else
162: if (p<0.0) SETERRQ(1,"IPNorm: The inner product is not well defined");
163: *norm = PetscSqrtScalar(p);
164: #endif
165: }
167: return(0);
168: }
172: /*@
173: IPInnerProduct - Computes the inner product of two vectors.
175: Collective on IP and Vec
177: Input Parameters:
178: + ip - the inner product context
179: . x - input vector
180: - y - input vector
182: Output Parameter:
183: . p - result of the inner product
185: Notes:
186: This function will usually compute the standard dot product of vectors
187: x and y, (x,y)=y^H x. However this behaviour may be different if changed
188: via IPSetBilinearForm(). This allows use of other inner products such as
189: the indefinite product y^T x for complex symmetric problems or the
190: B-inner product for positive definite B, (x,y)_B=y^H Bx.
192: Level: developer
194: .seealso: IPSetBilinearForm(), IPApplyB(), VecDot(), IPMInnerProduct()
195: @*/
196: PetscErrorCode IPInnerProduct(IP ip,Vec x,Vec y,PetscScalar *p)
197: {
205:
206: PetscLogEventBegin(IP_InnerProduct,ip,x,0,0);
207: ip->innerproducts++;
208: if (ip->matrix) {
209: IPApplyMatrix_Private(ip,x);
210: if (ip->bilinear_form == IPINNER_HERMITIAN) {
211: VecDot(ip->Bx,y,p);
212: } else {
213: VecTDot(ip->Bx,y,p);
214: }
215: } else {
216: if (ip->bilinear_form == IPINNER_HERMITIAN) {
217: VecDot(x,y,p);
218: } else {
219: VecTDot(x,y,p);
220: }
221: }
222: PetscLogEventEnd(IP_InnerProduct,ip,x,0,0);
223: return(0);
224: }
228: /*@
229: IPInnerProductBegin - Starts a split phase inner product computation.
231: Input Parameters:
232: + ip - the inner product context
233: . x - the first vector
234: . y - the second vector
235: - p - where the result will go
237: Level: developer
239: Notes:
240: Each call to IPInnerProductBegin() should be paired with a call to IPInnerProductEnd().
242: .seealso: IPInnerProductEnd(), IPInnerProduct(), IPNorm(), IPNormBegin(),
243: IPNormEnd(), IPMInnerProduct()
245: @*/
246: PetscErrorCode IPInnerProductBegin(IP ip,Vec x,Vec y,PetscScalar *p)
247: {
255:
256: PetscLogEventBegin(IP_InnerProduct,ip,x,0,0);
257: ip->innerproducts++;
258: if (ip->matrix) {
259: IPApplyMatrix_Private(ip,x);
260: if (ip->bilinear_form == IPINNER_HERMITIAN) {
261: VecDotBegin(ip->Bx,y,p);
262: } else {
263: VecTDotBegin(ip->Bx,y,p);
264: }
265: } else {
266: if (ip->bilinear_form == IPINNER_HERMITIAN) {
267: VecDotBegin(x,y,p);
268: } else {
269: VecTDotBegin(x,y,p);
270: }
271: }
272: PetscLogEventEnd(IP_InnerProduct,ip,x,0,0);
273: return(0);
274: }
278: /*@
279: IPInnerProductEnd - Ends a split phase inner product computation.
281: Input Parameters:
282: + ip - the inner product context
283: . x - the first vector
284: - y - the second vector
286: Output Parameter:
287: . p - result of the inner product
289: Level: developer
291: Notes:
292: Each call to IPInnerProductBegin() should be paired with a call to IPInnerProductEnd().
294: .seealso: IPInnerProductBegin(), IPInnerProduct(), IPNorm(), IPNormBegin(),
295: IPNormEnd(), IPMInnerProduct()
297: @*/
298: PetscErrorCode IPInnerProductEnd(IP ip,Vec x,Vec y,PetscScalar *p)
299: {
307:
308: PetscLogEventBegin(IP_InnerProduct,ip,x,0,0);
309: if (ip->matrix) {
310: if (ip->bilinear_form == IPINNER_HERMITIAN) {
311: VecDotEnd(ip->Bx,y,p);
312: } else {
313: VecTDotEnd(ip->Bx,y,p);
314: }
315: } else {
316: if (ip->bilinear_form == IPINNER_HERMITIAN) {
317: VecDotEnd(x,y,p);
318: } else {
319: VecTDotEnd(x,y,p);
320: }
321: }
322: PetscLogEventEnd(IP_InnerProduct,ip,x,0,0);
323: return(0);
324: }
328: /*@
329: IPMInnerProduct - Computes the inner products a vector x with a set of
330: vectors (columns of Y).
332: Collective on IP and Vec
334: Input Parameters:
335: + ip - the inner product context
336: . x - the first input vector
337: . n - number of vectors in y
338: - y - array of vectors
340: Output Parameter:
341: . p - result of the inner products
343: Notes:
344: This function will usually compute the standard dot product of x and y_i,
345: (x,y_i)=y_i^H x, for each column of Y. However this behaviour may be different
346: if changed via IPSetBilinearForm(). This allows use of other inner products
347: such as the indefinite product y_i^T x for complex symmetric problems or the
348: B-inner product for positive definite B, (x,y_i)_B=y_i^H Bx.
350: Level: developer
352: .seealso: IPSetBilinearForm(), IPApplyB(), VecMDot(), IPInnerProduct()
353: @*/
354: PetscErrorCode IPMInnerProduct(IP ip,Vec x,PetscInt n,const Vec y[],PetscScalar *p)
355: {
364:
365: PetscLogEventBegin(IP_InnerProduct,ip,x,0,0);
366: ip->innerproducts += n;
367: if (ip->matrix) {
368: IPApplyMatrix_Private(ip,x);
369: if (ip->bilinear_form == IPINNER_HERMITIAN) {
370: VecMDot(ip->Bx,n,y,p);
371: } else {
372: VecMTDot(ip->Bx,n,y,p);
373: }
374: } else {
375: if (ip->bilinear_form == IPINNER_HERMITIAN) {
376: VecMDot(x,n,y,p);
377: } else {
378: VecMTDot(x,n,y,p);
379: }
380: }
381: PetscLogEventEnd(IP_InnerProduct,ip,x,0,0);
382: return(0);
383: }
387: /*@
388: IPMInnerProductBegin - Starts a split phase multiple inner product computation.
390: Input Parameters:
391: + ip - the inner product context
392: . x - the first input vector
393: . n - number of vectors in y
394: . y - array of vectors
395: - p - where the result will go
397: Level: developer
399: Notes:
400: Each call to IPMInnerProductBegin() should be paired with a call to IPMInnerProductEnd().
402: .seealso: IPMInnerProductEnd(), IPMInnerProduct(), IPNorm(), IPNormBegin(),
403: IPNormEnd(), IPInnerProduct()
405: @*/
406: PetscErrorCode IPMInnerProductBegin(IP ip,Vec x,PetscInt n,const Vec y[],PetscScalar *p)
407: {
416:
417: PetscLogEventBegin(IP_InnerProduct,ip,x,0,0);
418: ip->innerproducts += n;
419: if (ip->matrix) {
420: IPApplyMatrix_Private(ip,x);
421: if (ip->bilinear_form == IPINNER_HERMITIAN) {
422: VecMDotBegin(ip->Bx,n,y,p);
423: } else {
424: VecMTDotBegin(ip->Bx,n,y,p);
425: }
426: } else {
427: if (ip->bilinear_form == IPINNER_HERMITIAN) {
428: VecMDotBegin(x,n,y,p);
429: } else {
430: VecMTDotBegin(x,n,y,p);
431: }
432: }
433: PetscLogEventEnd(IP_InnerProduct,ip,x,0,0);
434: return(0);
435: }
439: /*@
440: IPMInnerProductEnd - Ends a split phase multiple inner product computation.
442: Input Parameters:
443: + ip - the inner product context
444: . x - the first input vector
445: . n - number of vectors in y
446: - y - array of vectors
448: Output Parameter:
449: . p - result of the inner products
451: Level: developer
453: Notes:
454: Each call to IPMInnerProductBegin() should be paired with a call to IPMInnerProductEnd().
456: .seealso: IPMInnerProductBegin(), IPMInnerProduct(), IPNorm(), IPNormBegin(),
457: IPNormEnd(), IPInnerProduct()
459: @*/
460: PetscErrorCode IPMInnerProductEnd(IP ip,Vec x,PetscInt n,const Vec y[],PetscScalar *p)
461: {
470:
471: PetscLogEventBegin(IP_InnerProduct,ip,x,0,0);
472: if (ip->matrix) {
473: if (ip->bilinear_form == IPINNER_HERMITIAN) {
474: VecMDotEnd(ip->Bx,n,y,p);
475: } else {
476: VecMTDotEnd(ip->Bx,n,y,p);
477: }
478: } else {
479: if (ip->bilinear_form == IPINNER_HERMITIAN) {
480: VecMDotEnd(x,n,y,p);
481: } else {
482: VecMTDotEnd(x,n,y,p);
483: }
484: }
485: PetscLogEventEnd(IP_InnerProduct,ip,x,0,0);
486: return(0);
487: }