Actual source code: projection.c
1: #include <petsc/private/vecimpl.h>
3: /*@
4: VecWhichEqual - Creates an index set containing the indices
5: where the vectors Vec1 and Vec2 have identical elements.
7: Collective on Vec
9: Input Parameters:
10: . Vec1, Vec2 - the two vectors to compare
12: OutputParameter:
13: . S - The index set containing the indices i where vec1[i] == vec2[i]
15: Notes:
16: the two vectors must have the same parallel layout
18: Level: advanced
19: @*/
20: PetscErrorCode VecWhichEqual(Vec Vec1, Vec Vec2, IS *S)
21: {
22: PetscInt i, n_same = 0;
23: PetscInt n, low, high;
24: PetscInt *same = NULL;
25: const PetscScalar *v1, *v2;
27: PetscFunctionBegin;
30: PetscCheckSameComm(Vec1, 1, Vec2, 2);
31: VecCheckSameSize(Vec1, 1, Vec2, 2);
33: PetscCall(VecGetOwnershipRange(Vec1, &low, &high));
34: PetscCall(VecGetLocalSize(Vec1, &n));
35: if (n > 0) {
36: if (Vec1 == Vec2) {
37: PetscCall(VecGetArrayRead(Vec1, &v1));
38: v2 = v1;
39: } else {
40: PetscCall(VecGetArrayRead(Vec1, &v1));
41: PetscCall(VecGetArrayRead(Vec2, &v2));
42: }
44: PetscCall(PetscMalloc1(n, &same));
46: for (i = 0; i < n; ++i) {
47: if (v1[i] == v2[i]) {
48: same[n_same] = low + i;
49: ++n_same;
50: }
51: }
53: if (Vec1 == Vec2) {
54: PetscCall(VecRestoreArrayRead(Vec1, &v1));
55: } else {
56: PetscCall(VecRestoreArrayRead(Vec1, &v1));
57: PetscCall(VecRestoreArrayRead(Vec2, &v2));
58: }
59: }
60: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)Vec1), n_same, same, PETSC_OWN_POINTER, S));
61: PetscFunctionReturn(PETSC_SUCCESS);
62: }
64: /*@
65: VecWhichLessThan - Creates an index set containing the indices
66: where the vectors Vec1 < Vec2
68: Collective on S
70: Input Parameters:
71: . Vec1, Vec2 - the two vectors to compare
73: OutputParameter:
74: . S - The index set containing the indices i where vec1[i] < vec2[i]
76: Notes:
77: The two vectors must have the same parallel layout
79: For complex numbers this only compares the real part
81: Level: advanced
82: @*/
83: PetscErrorCode VecWhichLessThan(Vec Vec1, Vec Vec2, IS *S)
84: {
85: PetscInt i, n_lt = 0;
86: PetscInt n, low, high;
87: PetscInt *lt = NULL;
88: const PetscScalar *v1, *v2;
90: PetscFunctionBegin;
93: PetscCheckSameComm(Vec1, 1, Vec2, 2);
94: VecCheckSameSize(Vec1, 1, Vec2, 2);
96: PetscCall(VecGetOwnershipRange(Vec1, &low, &high));
97: PetscCall(VecGetLocalSize(Vec1, &n));
98: if (n > 0) {
99: if (Vec1 == Vec2) {
100: PetscCall(VecGetArrayRead(Vec1, &v1));
101: v2 = v1;
102: } else {
103: PetscCall(VecGetArrayRead(Vec1, &v1));
104: PetscCall(VecGetArrayRead(Vec2, &v2));
105: }
107: PetscCall(PetscMalloc1(n, <));
109: for (i = 0; i < n; ++i) {
110: if (PetscRealPart(v1[i]) < PetscRealPart(v2[i])) {
111: lt[n_lt] = low + i;
112: ++n_lt;
113: }
114: }
116: if (Vec1 == Vec2) {
117: PetscCall(VecRestoreArrayRead(Vec1, &v1));
118: } else {
119: PetscCall(VecRestoreArrayRead(Vec1, &v1));
120: PetscCall(VecRestoreArrayRead(Vec2, &v2));
121: }
122: }
123: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)Vec1), n_lt, lt, PETSC_OWN_POINTER, S));
124: PetscFunctionReturn(PETSC_SUCCESS);
125: }
127: /*@
128: VecWhichGreaterThan - Creates an index set containing the indices
129: where the vectors Vec1 > Vec2
131: Collective on S
133: Input Parameters:
134: . Vec1, Vec2 - the two vectors to compare
136: OutputParameter:
137: . S - The index set containing the indices i where vec1[i] > vec2[i]
139: Notes:
140: The two vectors must have the same parallel layout
142: For complex numbers this only compares the real part
144: Level: advanced
145: @*/
146: PetscErrorCode VecWhichGreaterThan(Vec Vec1, Vec Vec2, IS *S)
147: {
148: PetscInt i, n_gt = 0;
149: PetscInt n, low, high;
150: PetscInt *gt = NULL;
151: const PetscScalar *v1, *v2;
153: PetscFunctionBegin;
156: PetscCheckSameComm(Vec1, 1, Vec2, 2);
157: VecCheckSameSize(Vec1, 1, Vec2, 2);
159: PetscCall(VecGetOwnershipRange(Vec1, &low, &high));
160: PetscCall(VecGetLocalSize(Vec1, &n));
161: if (n > 0) {
162: if (Vec1 == Vec2) {
163: PetscCall(VecGetArrayRead(Vec1, &v1));
164: v2 = v1;
165: } else {
166: PetscCall(VecGetArrayRead(Vec1, &v1));
167: PetscCall(VecGetArrayRead(Vec2, &v2));
168: }
170: PetscCall(PetscMalloc1(n, >));
172: for (i = 0; i < n; ++i) {
173: if (PetscRealPart(v1[i]) > PetscRealPart(v2[i])) {
174: gt[n_gt] = low + i;
175: ++n_gt;
176: }
177: }
179: if (Vec1 == Vec2) {
180: PetscCall(VecRestoreArrayRead(Vec1, &v1));
181: } else {
182: PetscCall(VecRestoreArrayRead(Vec1, &v1));
183: PetscCall(VecRestoreArrayRead(Vec2, &v2));
184: }
185: }
186: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)Vec1), n_gt, gt, PETSC_OWN_POINTER, S));
187: PetscFunctionReturn(PETSC_SUCCESS);
188: }
190: /*@
191: VecWhichBetween - Creates an index set containing the indices
192: where VecLow < V < VecHigh
194: Collective on S
196: Input Parameters:
197: + VecLow - lower bound
198: . V - Vector to compare
199: - VecHigh - higher bound
201: OutputParameter:
202: . S - The index set containing the indices i where veclow[i] < v[i] < vechigh[i]
204: Notes:
205: The vectors must have the same parallel layout
207: For complex numbers this only compares the real part
209: Level: advanced
210: @*/
211: PetscErrorCode VecWhichBetween(Vec VecLow, Vec V, Vec VecHigh, IS *S)
212: {
213: PetscInt i, n_vm = 0;
214: PetscInt n, low, high;
215: PetscInt *vm = NULL;
216: const PetscScalar *v1, *v2, *vmiddle;
218: PetscFunctionBegin;
222: PetscCheckSameComm(V, 2, VecLow, 1);
223: PetscCheckSameComm(V, 2, VecHigh, 3);
224: VecCheckSameSize(V, 2, VecLow, 1);
225: VecCheckSameSize(V, 2, VecHigh, 3);
227: PetscCall(VecGetOwnershipRange(VecLow, &low, &high));
228: PetscCall(VecGetLocalSize(VecLow, &n));
229: if (n > 0) {
230: PetscCall(VecGetArrayRead(VecLow, &v1));
231: if (VecLow != VecHigh) {
232: PetscCall(VecGetArrayRead(VecHigh, &v2));
233: } else {
234: v2 = v1;
235: }
236: if (V != VecLow && V != VecHigh) {
237: PetscCall(VecGetArrayRead(V, &vmiddle));
238: } else if (V == VecLow) {
239: vmiddle = v1;
240: } else {
241: vmiddle = v2;
242: }
244: PetscCall(PetscMalloc1(n, &vm));
246: for (i = 0; i < n; ++i) {
247: if (PetscRealPart(v1[i]) < PetscRealPart(vmiddle[i]) && PetscRealPart(vmiddle[i]) < PetscRealPart(v2[i])) {
248: vm[n_vm] = low + i;
249: ++n_vm;
250: }
251: }
253: PetscCall(VecRestoreArrayRead(VecLow, &v1));
254: if (VecLow != VecHigh) PetscCall(VecRestoreArrayRead(VecHigh, &v2));
255: if (V != VecLow && V != VecHigh) PetscCall(VecRestoreArrayRead(V, &vmiddle));
256: }
257: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)V), n_vm, vm, PETSC_OWN_POINTER, S));
258: PetscFunctionReturn(PETSC_SUCCESS);
259: }
261: /*@
262: VecWhichBetweenOrEqual - Creates an index set containing the indices
263: where VecLow <= V <= VecHigh
265: Collective on S
267: Input Parameters:
268: + VecLow - lower bound
269: . V - Vector to compare
270: - VecHigh - higher bound
272: OutputParameter:
273: . S - The index set containing the indices i where veclow[i] <= v[i] <= vechigh[i]
275: Level: advanced
276: @*/
278: PetscErrorCode VecWhichBetweenOrEqual(Vec VecLow, Vec V, Vec VecHigh, IS *S)
279: {
280: PetscInt i, n_vm = 0;
281: PetscInt n, low, high;
282: PetscInt *vm = NULL;
283: const PetscScalar *v1, *v2, *vmiddle;
285: PetscFunctionBegin;
289: PetscCheckSameComm(V, 2, VecLow, 1);
290: PetscCheckSameComm(V, 2, VecHigh, 3);
291: VecCheckSameSize(V, 2, VecLow, 1);
292: VecCheckSameSize(V, 2, VecHigh, 3);
294: PetscCall(VecGetOwnershipRange(VecLow, &low, &high));
295: PetscCall(VecGetLocalSize(VecLow, &n));
296: if (n > 0) {
297: PetscCall(VecGetArrayRead(VecLow, &v1));
298: if (VecLow != VecHigh) {
299: PetscCall(VecGetArrayRead(VecHigh, &v2));
300: } else {
301: v2 = v1;
302: }
303: if (V != VecLow && V != VecHigh) {
304: PetscCall(VecGetArrayRead(V, &vmiddle));
305: } else if (V == VecLow) {
306: vmiddle = v1;
307: } else {
308: vmiddle = v2;
309: }
311: PetscCall(PetscMalloc1(n, &vm));
313: for (i = 0; i < n; ++i) {
314: if (PetscRealPart(v1[i]) <= PetscRealPart(vmiddle[i]) && PetscRealPart(vmiddle[i]) <= PetscRealPart(v2[i])) {
315: vm[n_vm] = low + i;
316: ++n_vm;
317: }
318: }
320: PetscCall(VecRestoreArrayRead(VecLow, &v1));
321: if (VecLow != VecHigh) PetscCall(VecRestoreArrayRead(VecHigh, &v2));
322: if (V != VecLow && V != VecHigh) PetscCall(VecRestoreArrayRead(V, &vmiddle));
323: }
324: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)V), n_vm, vm, PETSC_OWN_POINTER, S));
325: PetscFunctionReturn(PETSC_SUCCESS);
326: }
328: /*@
329: VecWhichInactive - Creates an index set containing the indices
330: where one of the following holds:
331: a) VecLow(i) < V(i) < VecHigh(i)
332: b) VecLow(i) = V(i) and D(i) <= 0 (< 0 when Strong is true)
333: c) VecHigh(i) = V(i) and D(i) >= 0 (> 0 when Strong is true)
335: Collective on S
337: Input Parameters:
338: + VecLow - lower bound
339: . V - Vector to compare
340: . D - Direction to compare
341: . VecHigh - higher bound
342: - Strong - indicator for applying strongly inactive test
344: OutputParameter:
345: . S - The index set containing the indices i where the bound is inactive
347: Level: advanced
348: @*/
350: PetscErrorCode VecWhichInactive(Vec VecLow, Vec V, Vec D, Vec VecHigh, PetscBool Strong, IS *S)
351: {
352: PetscInt i, n_vm = 0;
353: PetscInt n, low, high;
354: PetscInt *vm = NULL;
355: const PetscScalar *v1, *v2, *v, *d;
357: PetscFunctionBegin;
358: if (!VecLow && !VecHigh) {
359: PetscCall(VecGetOwnershipRange(V, &low, &high));
360: PetscCall(ISCreateStride(PetscObjectComm((PetscObject)V), high - low, low, 1, S));
361: PetscFunctionReturn(PETSC_SUCCESS);
362: }
367: PetscCheckSameComm(V, 2, VecLow, 1);
368: PetscCheckSameComm(V, 2, D, 3);
369: PetscCheckSameComm(V, 2, VecHigh, 4);
370: VecCheckSameSize(V, 2, VecLow, 1);
371: VecCheckSameSize(V, 2, D, 3);
372: VecCheckSameSize(V, 2, VecHigh, 4);
374: PetscCall(VecGetOwnershipRange(VecLow, &low, &high));
375: PetscCall(VecGetLocalSize(VecLow, &n));
376: if (n > 0) {
377: PetscCall(VecGetArrayRead(VecLow, &v1));
378: if (VecLow != VecHigh) {
379: PetscCall(VecGetArrayRead(VecHigh, &v2));
380: } else {
381: v2 = v1;
382: }
383: if (V != VecLow && V != VecHigh) {
384: PetscCall(VecGetArrayRead(V, &v));
385: } else if (V == VecLow) {
386: v = v1;
387: } else {
388: v = v2;
389: }
390: if (D != VecLow && D != VecHigh && D != V) {
391: PetscCall(VecGetArrayRead(D, &d));
392: } else if (D == VecLow) {
393: d = v1;
394: } else if (D == VecHigh) {
395: d = v2;
396: } else {
397: d = v;
398: }
400: PetscCall(PetscMalloc1(n, &vm));
402: if (Strong) {
403: for (i = 0; i < n; ++i) {
404: if (PetscRealPart(v1[i]) < PetscRealPart(v[i]) && PetscRealPart(v[i]) < PetscRealPart(v2[i])) {
405: vm[n_vm] = low + i;
406: ++n_vm;
407: } else if (PetscRealPart(v1[i]) == PetscRealPart(v[i]) && PetscRealPart(d[i]) < 0) {
408: vm[n_vm] = low + i;
409: ++n_vm;
410: } else if (PetscRealPart(v2[i]) == PetscRealPart(v[i]) && PetscRealPart(d[i]) > 0) {
411: vm[n_vm] = low + i;
412: ++n_vm;
413: }
414: }
415: } else {
416: for (i = 0; i < n; ++i) {
417: if (PetscRealPart(v1[i]) < PetscRealPart(v[i]) && PetscRealPart(v[i]) < PetscRealPart(v2[i])) {
418: vm[n_vm] = low + i;
419: ++n_vm;
420: } else if (PetscRealPart(v1[i]) == PetscRealPart(v[i]) && PetscRealPart(d[i]) <= 0) {
421: vm[n_vm] = low + i;
422: ++n_vm;
423: } else if (PetscRealPart(v2[i]) == PetscRealPart(v[i]) && PetscRealPart(d[i]) >= 0) {
424: vm[n_vm] = low + i;
425: ++n_vm;
426: }
427: }
428: }
430: PetscCall(VecRestoreArrayRead(VecLow, &v1));
431: if (VecLow != VecHigh) PetscCall(VecRestoreArrayRead(VecHigh, &v2));
432: if (V != VecLow && V != VecHigh) PetscCall(VecRestoreArrayRead(V, &v));
433: if (D != VecLow && D != VecHigh && D != V) PetscCall(VecRestoreArrayRead(D, &d));
434: }
435: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)V), n_vm, vm, PETSC_OWN_POINTER, S));
436: PetscFunctionReturn(PETSC_SUCCESS);
437: }
439: /*@
440: VecISAXPY - Adds a reduced vector to the appropriate elements of a full-space vector.
441: vfull[is[i]] += alpha*vreduced[i]
443: Input Parameters:
444: + vfull - the full-space vector
445: . is - the index set for the reduced space
446: . alpha - the scalar coefficient
447: - vreduced - the reduced-space vector
449: Output Parameters:
450: . vfull - the sum of the full-space vector and reduced-space vector
452: Notes:
453: The index set identifies entries in the global vector.
454: Negative indices are skipped; indices outside the ownership range of vfull will raise an error.
456: Level: advanced
458: .seealso: `VecISCopy()`, `VecISSet()`, `VecAXPY()`
459: @*/
460: PetscErrorCode VecISAXPY(Vec vfull, IS is, PetscScalar alpha, Vec vreduced)
461: {
462: PetscInt nfull, nreduced;
464: PetscFunctionBegin;
468: PetscCall(VecGetSize(vfull, &nfull));
469: PetscCall(VecGetSize(vreduced, &nreduced));
471: if (nfull == nreduced) { /* Also takes care of masked vectors */
472: PetscCall(VecAXPY(vfull, alpha, vreduced));
473: } else {
474: PetscScalar *y;
475: const PetscScalar *x;
476: PetscInt i, n, m, rstart, rend;
477: const PetscInt *id;
479: PetscCall(VecGetArray(vfull, &y));
480: PetscCall(VecGetArrayRead(vreduced, &x));
481: PetscCall(ISGetIndices(is, &id));
482: PetscCall(ISGetLocalSize(is, &n));
483: PetscCall(VecGetLocalSize(vreduced, &m));
484: PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_SUP, "IS local length not equal to Vec local length");
485: PetscCall(VecGetOwnershipRange(vfull, &rstart, &rend));
486: y -= rstart;
487: if (alpha == 1.0) {
488: for (i = 0; i < n; ++i) {
489: if (id[i] < 0) continue;
490: PetscCheck(id[i] >= rstart && id[i] < rend, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only owned values supported");
491: y[id[i]] += x[i];
492: }
493: } else {
494: for (i = 0; i < n; ++i) {
495: if (id[i] < 0) continue;
496: PetscCheck(id[i] >= rstart && id[i] < rend, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only owned values supported");
497: y[id[i]] += alpha * x[i];
498: }
499: }
500: y += rstart;
501: PetscCall(ISRestoreIndices(is, &id));
502: PetscCall(VecRestoreArray(vfull, &y));
503: PetscCall(VecRestoreArrayRead(vreduced, &x));
504: }
505: PetscFunctionReturn(PETSC_SUCCESS);
506: }
508: /*@
509: VecISCopy - Copies between a reduced vector and the appropriate elements of a full-space vector.
511: Input Parameters:
512: + vfull - the full-space vector
513: . is - the index set for the reduced space
514: . mode - the direction of copying, SCATTER_FORWARD or SCATTER_REVERSE
515: - vreduced - the reduced-space vector
517: Output Parameters:
518: . vfull - the sum of the full-space vector and reduced-space vector
520: Notes:
521: The index set identifies entries in the global vector.
522: Negative indices are skipped; indices outside the ownership range of vfull will raise an error.
524: mode == SCATTER_FORWARD: vfull[is[i]] = vreduced[i]
525: mode == SCATTER_REVERSE: vreduced[i] = vfull[is[i]]
527: Level: advanced
529: .seealso: `VecISSet()`, `VecISAXPY()`, `VecCopy()`
530: @*/
531: PetscErrorCode VecISCopy(Vec vfull, IS is, ScatterMode mode, Vec vreduced)
532: {
533: PetscInt nfull, nreduced;
535: PetscFunctionBegin;
539: PetscCall(VecGetSize(vfull, &nfull));
540: PetscCall(VecGetSize(vreduced, &nreduced));
542: if (nfull == nreduced) { /* Also takes care of masked vectors */
543: if (mode == SCATTER_FORWARD) {
544: PetscCall(VecCopy(vreduced, vfull));
545: } else {
546: PetscCall(VecCopy(vfull, vreduced));
547: }
548: } else {
549: const PetscInt *id;
550: PetscInt i, n, m, rstart, rend;
552: PetscCall(ISGetIndices(is, &id));
553: PetscCall(ISGetLocalSize(is, &n));
554: PetscCall(VecGetLocalSize(vreduced, &m));
555: PetscCall(VecGetOwnershipRange(vfull, &rstart, &rend));
556: PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_SUP, "IS local length %" PetscInt_FMT " not equal to Vec local length %" PetscInt_FMT, n, m);
557: if (mode == SCATTER_FORWARD) {
558: PetscScalar *y;
559: const PetscScalar *x;
561: PetscCall(VecGetArray(vfull, &y));
562: PetscCall(VecGetArrayRead(vreduced, &x));
563: y -= rstart;
564: for (i = 0; i < n; ++i) {
565: if (id[i] < 0) continue;
566: PetscCheck(id[i] >= rstart && id[i] < rend, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only owned values supported");
567: y[id[i]] = x[i];
568: }
569: y += rstart;
570: PetscCall(VecRestoreArrayRead(vreduced, &x));
571: PetscCall(VecRestoreArray(vfull, &y));
572: } else if (mode == SCATTER_REVERSE) {
573: PetscScalar *x;
574: const PetscScalar *y;
576: PetscCall(VecGetArrayRead(vfull, &y));
577: PetscCall(VecGetArray(vreduced, &x));
578: for (i = 0; i < n; ++i) {
579: if (id[i] < 0) continue;
580: PetscCheck(id[i] >= rstart && id[i] < rend, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only owned values supported");
581: x[i] = y[id[i] - rstart];
582: }
583: PetscCall(VecRestoreArray(vreduced, &x));
584: PetscCall(VecRestoreArrayRead(vfull, &y));
585: } else SETERRQ(PetscObjectComm((PetscObject)vfull), PETSC_ERR_ARG_WRONG, "Only forward or reverse modes are legal");
586: PetscCall(ISRestoreIndices(is, &id));
587: }
588: PetscFunctionReturn(PETSC_SUCCESS);
589: }
591: /*@
592: ISComplementVec - Creates the complement of the index set relative to a layout defined by a Vec
594: Collective on IS
596: Input Parameters:
597: + S - a PETSc IS
598: - V - the reference vector space
600: Output Parameter:
601: . T - the complement of S
603: Level: advanced
605: .seealso: `ISCreateGeneral()`
606: @*/
607: PetscErrorCode ISComplementVec(IS S, Vec V, IS *T)
608: {
609: PetscInt start, end;
611: PetscFunctionBegin;
612: PetscCall(VecGetOwnershipRange(V, &start, &end));
613: PetscCall(ISComplement(S, start, end, T));
614: PetscFunctionReturn(PETSC_SUCCESS);
615: }
617: /*@
618: VecISSet - Sets the elements of a vector, specified by an index set, to a constant
620: Input Parameters:
621: + V - the vector
622: . S - index set for the locations in the vector
623: - c - the constant
625: Notes:
626: The index set identifies entries in the global vector.
627: Negative indices are skipped; indices outside the ownership range of V will raise an error.
629: Level: advanced
631: .seealso: `VecISCopy()`, `VecISAXPY()`, `VecSet()`
632: @*/
633: PetscErrorCode VecISSet(Vec V, IS S, PetscScalar c)
634: {
635: PetscInt nloc, low, high, i;
636: const PetscInt *s;
637: PetscScalar *v;
639: PetscFunctionBegin;
640: if (!S) PetscFunctionReturn(PETSC_SUCCESS); /* simply return with no-op if the index set is NULL */
645: PetscCall(VecGetOwnershipRange(V, &low, &high));
646: PetscCall(ISGetLocalSize(S, &nloc));
647: PetscCall(ISGetIndices(S, &s));
648: PetscCall(VecGetArray(V, &v));
649: for (i = 0; i < nloc; ++i) {
650: if (s[i] < 0) continue;
651: PetscCheck(s[i] >= low && s[i] < high, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only owned values supported");
652: v[s[i] - low] = c;
653: }
654: PetscCall(ISRestoreIndices(S, &s));
655: PetscCall(VecRestoreArray(V, &v));
656: PetscFunctionReturn(PETSC_SUCCESS);
657: }
659: #if !defined(PETSC_USE_COMPLEX)
660: /*@C
661: VecBoundGradientProjection - Projects vector according to this definition.
662: If XL[i] < X[i] < XU[i], then GP[i] = G[i];
663: If X[i] <= XL[i], then GP[i] = min(G[i],0);
664: If X[i] >= XU[i], then GP[i] = max(G[i],0);
666: Input Parameters:
667: + G - current gradient vector
668: . X - current solution vector with XL[i] <= X[i] <= XU[i]
669: . XL - lower bounds
670: - XU - upper bounds
672: Output Parameter:
673: . GP - gradient projection vector
675: Notes:
676: GP may be the same vector as G
678: Level: advanced
679: @*/
680: PetscErrorCode VecBoundGradientProjection(Vec G, Vec X, Vec XL, Vec XU, Vec GP)
681: {
682: PetscInt n, i;
683: const PetscReal *xptr, *xlptr, *xuptr;
684: PetscReal *gptr, *gpptr;
685: PetscReal xval, gpval;
687: /* Project variables at the lower and upper bound */
688: PetscFunctionBegin;
694: if (!XL && !XU) {
695: PetscCall(VecCopy(G, GP));
696: PetscFunctionReturn(PETSC_SUCCESS);
697: }
699: PetscCall(VecGetLocalSize(X, &n));
701: PetscCall(VecGetArrayRead(X, &xptr));
702: PetscCall(VecGetArrayRead(XL, &xlptr));
703: PetscCall(VecGetArrayRead(XU, &xuptr));
704: PetscCall(VecGetArrayPair(G, GP, &gptr, &gpptr));
706: for (i = 0; i < n; ++i) {
707: gpval = gptr[i];
708: xval = xptr[i];
709: if (gpval > 0.0 && xval <= xlptr[i]) {
710: gpval = 0.0;
711: } else if (gpval < 0.0 && xval >= xuptr[i]) {
712: gpval = 0.0;
713: }
714: gpptr[i] = gpval;
715: }
717: PetscCall(VecRestoreArrayRead(X, &xptr));
718: PetscCall(VecRestoreArrayRead(XL, &xlptr));
719: PetscCall(VecRestoreArrayRead(XU, &xuptr));
720: PetscCall(VecRestoreArrayPair(G, GP, &gptr, &gpptr));
721: PetscFunctionReturn(PETSC_SUCCESS);
722: }
723: #endif
725: /*@
726: VecStepMaxBounded - See below
728: Collective on Vec
730: Input Parameters:
731: + X - vector with no negative entries
732: . XL - lower bounds
733: . XU - upper bounds
734: - DX - step direction, can have negative, positive or zero entries
736: Output Parameter:
737: . stepmax - minimum value so that X[i] + stepmax*DX[i] <= XL[i] or XU[i] <= X[i] + stepmax*DX[i]
739: Level: intermediate
741: @*/
742: PetscErrorCode VecStepMaxBounded(Vec X, Vec DX, Vec XL, Vec XU, PetscReal *stepmax)
743: {
744: PetscInt i, nn;
745: const PetscScalar *xx, *dx, *xl, *xu;
746: PetscReal localmax = 0;
748: PetscFunctionBegin;
754: PetscCall(VecGetArrayRead(X, &xx));
755: PetscCall(VecGetArrayRead(XL, &xl));
756: PetscCall(VecGetArrayRead(XU, &xu));
757: PetscCall(VecGetArrayRead(DX, &dx));
758: PetscCall(VecGetLocalSize(X, &nn));
759: for (i = 0; i < nn; i++) {
760: if (PetscRealPart(dx[i]) > 0) {
761: localmax = PetscMax(localmax, PetscRealPart((xu[i] - xx[i]) / dx[i]));
762: } else if (PetscRealPart(dx[i]) < 0) {
763: localmax = PetscMax(localmax, PetscRealPart((xl[i] - xx[i]) / dx[i]));
764: }
765: }
766: PetscCall(VecRestoreArrayRead(X, &xx));
767: PetscCall(VecRestoreArrayRead(XL, &xl));
768: PetscCall(VecRestoreArrayRead(XU, &xu));
769: PetscCall(VecRestoreArrayRead(DX, &dx));
770: PetscCall(MPIU_Allreduce(&localmax, stepmax, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)X)));
771: PetscFunctionReturn(PETSC_SUCCESS);
772: }
774: /*@
775: VecStepBoundInfo - See below
777: Collective on Vec
779: Input Parameters:
780: + X - vector with no negative entries
781: . XL - lower bounds
782: . XU - upper bounds
783: - DX - step direction, can have negative, positive or zero entries
785: Output Parameters:
786: + boundmin - (may be NULL this it is not computed) maximum value so that XL[i] <= X[i] + boundmax*DX[i] <= XU[i]
787: . wolfemin - (may be NULL this it is not computed)
788: - boundmax - (may be NULL this it is not computed) minimum value so that X[i] + boundmax*DX[i] <= XL[i] or XU[i] <= X[i] + boundmax*DX[i]
790: Notes:
791: For complex numbers only compares the real part
793: Level: advanced
794: @*/
795: PetscErrorCode VecStepBoundInfo(Vec X, Vec DX, Vec XL, Vec XU, PetscReal *boundmin, PetscReal *wolfemin, PetscReal *boundmax)
796: {
797: PetscInt n, i;
798: const PetscScalar *x, *xl, *xu, *dx;
799: PetscReal t;
800: PetscReal localmin = PETSC_INFINITY, localwolfemin = PETSC_INFINITY, localmax = -1;
801: MPI_Comm comm;
803: PetscFunctionBegin;
809: PetscCall(VecGetArrayRead(X, &x));
810: PetscCall(VecGetArrayRead(XL, &xl));
811: PetscCall(VecGetArrayRead(XU, &xu));
812: PetscCall(VecGetArrayRead(DX, &dx));
813: PetscCall(VecGetLocalSize(X, &n));
814: for (i = 0; i < n; ++i) {
815: if (PetscRealPart(dx[i]) > 0 && PetscRealPart(xu[i]) < PETSC_INFINITY) {
816: t = PetscRealPart((xu[i] - x[i]) / dx[i]);
817: localmin = PetscMin(t, localmin);
818: if (localmin > 0) localwolfemin = PetscMin(t, localwolfemin);
819: localmax = PetscMax(t, localmax);
820: } else if (PetscRealPart(dx[i]) < 0 && PetscRealPart(xl[i]) > PETSC_NINFINITY) {
821: t = PetscRealPart((xl[i] - x[i]) / dx[i]);
822: localmin = PetscMin(t, localmin);
823: if (localmin > 0) localwolfemin = PetscMin(t, localwolfemin);
824: localmax = PetscMax(t, localmax);
825: }
826: }
828: PetscCall(VecRestoreArrayRead(X, &x));
829: PetscCall(VecRestoreArrayRead(XL, &xl));
830: PetscCall(VecRestoreArrayRead(XU, &xu));
831: PetscCall(VecRestoreArrayRead(DX, &dx));
832: PetscCall(PetscObjectGetComm((PetscObject)X, &comm));
834: if (boundmin) {
835: PetscCall(MPIU_Allreduce(&localmin, boundmin, 1, MPIU_REAL, MPIU_MIN, comm));
836: PetscCall(PetscInfo(X, "Step Bound Info: Closest Bound: %20.19e\n", (double)*boundmin));
837: }
838: if (wolfemin) {
839: PetscCall(MPIU_Allreduce(&localwolfemin, wolfemin, 1, MPIU_REAL, MPIU_MIN, comm));
840: PetscCall(PetscInfo(X, "Step Bound Info: Wolfe: %20.19e\n", (double)*wolfemin));
841: }
842: if (boundmax) {
843: PetscCall(MPIU_Allreduce(&localmax, boundmax, 1, MPIU_REAL, MPIU_MAX, comm));
844: if (*boundmax < 0) *boundmax = PETSC_INFINITY;
845: PetscCall(PetscInfo(X, "Step Bound Info: Max: %20.19e\n", (double)*boundmax));
846: }
847: PetscFunctionReturn(PETSC_SUCCESS);
848: }
850: /*@
851: VecStepMax - Returns the largest value so that x[i] + step*DX[i] >= 0 for all i
853: Collective on Vec
855: Input Parameters:
856: + X - vector with no negative entries
857: - DX - a step direction, can have negative, positive or zero entries
859: Output Parameter:
860: . step - largest value such that x[i] + step*DX[i] >= 0 for all i
862: Notes:
863: For complex numbers only compares the real part
865: Level: advanced
866: @*/
867: PetscErrorCode VecStepMax(Vec X, Vec DX, PetscReal *step)
868: {
869: PetscInt i, nn;
870: PetscReal stepmax = PETSC_INFINITY;
871: const PetscScalar *xx, *dx;
873: PetscFunctionBegin;
877: PetscCall(VecGetLocalSize(X, &nn));
878: PetscCall(VecGetArrayRead(X, &xx));
879: PetscCall(VecGetArrayRead(DX, &dx));
880: for (i = 0; i < nn; ++i) {
881: PetscCheck(PetscRealPart(xx[i]) >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vector must be positive");
882: if (PetscRealPart(dx[i]) < 0) stepmax = PetscMin(stepmax, PetscRealPart(-xx[i] / dx[i]));
883: }
884: PetscCall(VecRestoreArrayRead(X, &xx));
885: PetscCall(VecRestoreArrayRead(DX, &dx));
886: PetscCall(MPIU_Allreduce(&stepmax, step, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)X)));
887: PetscFunctionReturn(PETSC_SUCCESS);
888: }
890: /*@
891: VecPow - Replaces each component of a vector by x_i^p
893: Logically Collective on v
895: Input Parameters:
896: + v - the vector
897: - p - the exponent to use on each element
899: Level: intermediate
901: @*/
902: PetscErrorCode VecPow(Vec v, PetscScalar p)
903: {
904: PetscInt n, i;
905: PetscScalar *v1;
907: PetscFunctionBegin;
910: PetscCall(VecGetArray(v, &v1));
911: PetscCall(VecGetLocalSize(v, &n));
913: if (1.0 == p) {
914: } else if (-1.0 == p) {
915: for (i = 0; i < n; ++i) v1[i] = 1.0 / v1[i];
916: } else if (0.0 == p) {
917: for (i = 0; i < n; ++i) {
918: /* Not-a-number left alone
919: Infinity set to one */
920: if (v1[i] == v1[i]) v1[i] = 1.0;
921: }
922: } else if (0.5 == p) {
923: for (i = 0; i < n; ++i) {
924: if (PetscRealPart(v1[i]) >= 0) {
925: v1[i] = PetscSqrtScalar(v1[i]);
926: } else {
927: v1[i] = PETSC_INFINITY;
928: }
929: }
930: } else if (-0.5 == p) {
931: for (i = 0; i < n; ++i) {
932: if (PetscRealPart(v1[i]) >= 0) {
933: v1[i] = 1.0 / PetscSqrtScalar(v1[i]);
934: } else {
935: v1[i] = PETSC_INFINITY;
936: }
937: }
938: } else if (2.0 == p) {
939: for (i = 0; i < n; ++i) v1[i] *= v1[i];
940: } else if (-2.0 == p) {
941: for (i = 0; i < n; ++i) v1[i] = 1.0 / (v1[i] * v1[i]);
942: } else {
943: for (i = 0; i < n; ++i) {
944: if (PetscRealPart(v1[i]) >= 0) {
945: v1[i] = PetscPowScalar(v1[i], p);
946: } else {
947: v1[i] = PETSC_INFINITY;
948: }
949: }
950: }
951: PetscCall(VecRestoreArray(v, &v1));
952: PetscFunctionReturn(PETSC_SUCCESS);
953: }
955: /*@
956: VecMedian - Computes the componentwise median of three vectors
957: and stores the result in this vector. Used primarily for projecting
958: a vector within upper and lower bounds.
960: Logically Collective
962: Input Parameters:
963: + Vec1 - The first vector
964: . Vec2 - The second vector
965: - Vec3 - The third vector
967: Output Parameter:
968: . VMedian - The median vector (this can be any one of the input vectors)
970: Level: advanced
971: @*/
972: PetscErrorCode VecMedian(Vec Vec1, Vec Vec2, Vec Vec3, Vec VMedian)
973: {
974: PetscInt i, n, low1, high1;
975: const PetscScalar *v1, *v2, *v3;
976: PetscScalar *vmed;
978: PetscFunctionBegin;
984: if (!Vec1 && !Vec3) {
985: PetscCall(VecCopy(Vec2, VMedian));
986: PetscFunctionReturn(PETSC_SUCCESS);
987: }
988: if (Vec1 == Vec2 || Vec1 == Vec3) {
989: PetscCall(VecCopy(Vec1, VMedian));
990: PetscFunctionReturn(PETSC_SUCCESS);
991: }
992: if (Vec2 == Vec3) {
993: PetscCall(VecCopy(Vec2, VMedian));
994: PetscFunctionReturn(PETSC_SUCCESS);
995: }
997: /* Assert that Vec1 != Vec2 and Vec2 != Vec3 */
1002: PetscCheckSameType(Vec1, 1, Vec2, 2);
1003: PetscCheckSameType(Vec1, 1, Vec3, 3);
1004: PetscCheckSameType(Vec1, 1, VMedian, 4);
1005: PetscCheckSameComm(Vec1, 1, Vec2, 2);
1006: PetscCheckSameComm(Vec1, 1, Vec3, 3);
1007: PetscCheckSameComm(Vec1, 1, VMedian, 4);
1008: VecCheckSameSize(Vec1, 1, Vec2, 2);
1009: VecCheckSameSize(Vec1, 1, Vec3, 3);
1010: VecCheckSameSize(Vec1, 1, VMedian, 4);
1012: PetscCall(VecGetOwnershipRange(Vec1, &low1, &high1));
1013: PetscCall(VecGetLocalSize(Vec1, &n));
1014: if (n > 0) {
1015: PetscCall(VecGetArray(VMedian, &vmed));
1016: if (Vec1 != VMedian) {
1017: PetscCall(VecGetArrayRead(Vec1, &v1));
1018: } else {
1019: v1 = vmed;
1020: }
1021: if (Vec2 != VMedian) {
1022: PetscCall(VecGetArrayRead(Vec2, &v2));
1023: } else {
1024: v2 = vmed;
1025: }
1026: if (Vec3 != VMedian) {
1027: PetscCall(VecGetArrayRead(Vec3, &v3));
1028: } else {
1029: v3 = vmed;
1030: }
1032: for (i = 0; i < n; ++i) vmed[i] = PetscMax(PetscMax(PetscMin(PetscRealPart(v1[i]), PetscRealPart(v2[i])), PetscMin(PetscRealPart(v1[i]), PetscRealPart(v3[i]))), PetscMin(PetscRealPart(v2[i]), PetscRealPart(v3[i])));
1034: PetscCall(VecRestoreArray(VMedian, &vmed));
1035: if (VMedian != Vec1) PetscCall(VecRestoreArrayRead(Vec1, &v1));
1036: if (VMedian != Vec2) PetscCall(VecRestoreArrayRead(Vec2, &v2));
1037: if (VMedian != Vec3) PetscCall(VecRestoreArrayRead(Vec3, &v3));
1038: }
1039: PetscFunctionReturn(PETSC_SUCCESS);
1040: }