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
9: Input Parameters:
10: + Vec1 - the first vector to compare
11: - Vec2 - the second two vector to compare
13: Output Parameter:
14: . S - The index set containing the indices i where vec1[i] == vec2[i]
16: Level: advanced
18: Note:
19: The two vectors must have the same parallel layout
21: .seealso: `Vec`
22: @*/
23: PetscErrorCode VecWhichEqual(Vec Vec1, Vec Vec2, IS *S)
24: {
25: PetscInt i, n_same = 0;
26: PetscInt n, low, high;
27: PetscInt *same = NULL;
28: const PetscScalar *v1, *v2;
30: PetscFunctionBegin;
33: PetscCheckSameComm(Vec1, 1, Vec2, 2);
34: VecCheckSameSize(Vec1, 1, Vec2, 2);
36: PetscCall(VecGetOwnershipRange(Vec1, &low, &high));
37: PetscCall(VecGetLocalSize(Vec1, &n));
38: if (n > 0) {
39: if (Vec1 == Vec2) {
40: PetscCall(VecGetArrayRead(Vec1, &v1));
41: v2 = v1;
42: } else {
43: PetscCall(VecGetArrayRead(Vec1, &v1));
44: PetscCall(VecGetArrayRead(Vec2, &v2));
45: }
47: PetscCall(PetscMalloc1(n, &same));
49: for (i = 0; i < n; ++i) {
50: if (v1[i] == v2[i]) {
51: same[n_same] = low + i;
52: ++n_same;
53: }
54: }
56: if (Vec1 == Vec2) {
57: PetscCall(VecRestoreArrayRead(Vec1, &v1));
58: } else {
59: PetscCall(VecRestoreArrayRead(Vec1, &v1));
60: PetscCall(VecRestoreArrayRead(Vec2, &v2));
61: }
62: }
63: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)Vec1), n_same, same, PETSC_OWN_POINTER, S));
64: PetscFunctionReturn(PETSC_SUCCESS);
65: }
67: /*@
68: VecWhichLessThan - Creates an index set containing the indices
69: where the vectors `Vec1` < `Vec2`
71: Collective
73: Input Parameters:
74: + Vec1 - the first vector to compare
75: - Vec2 - the second vector to compare
77: Output Parameter:
78: . S - The index set containing the indices i where vec1[i] < vec2[i]
80: Level: advanced
82: Notes:
83: The two vectors must have the same parallel layout
85: For complex numbers this only compares the real part
87: .seealso: `Vec`
88: @*/
89: PetscErrorCode VecWhichLessThan(Vec Vec1, Vec Vec2, IS *S)
90: {
91: PetscInt i, n_lt = 0;
92: PetscInt n, low, high;
93: PetscInt *lt = NULL;
94: const PetscScalar *v1, *v2;
96: PetscFunctionBegin;
99: PetscCheckSameComm(Vec1, 1, Vec2, 2);
100: VecCheckSameSize(Vec1, 1, Vec2, 2);
102: PetscCall(VecGetOwnershipRange(Vec1, &low, &high));
103: PetscCall(VecGetLocalSize(Vec1, &n));
104: if (n > 0) {
105: if (Vec1 == Vec2) {
106: PetscCall(VecGetArrayRead(Vec1, &v1));
107: v2 = v1;
108: } else {
109: PetscCall(VecGetArrayRead(Vec1, &v1));
110: PetscCall(VecGetArrayRead(Vec2, &v2));
111: }
113: PetscCall(PetscMalloc1(n, <));
115: for (i = 0; i < n; ++i) {
116: if (PetscRealPart(v1[i]) < PetscRealPart(v2[i])) {
117: lt[n_lt] = low + i;
118: ++n_lt;
119: }
120: }
122: if (Vec1 == Vec2) {
123: PetscCall(VecRestoreArrayRead(Vec1, &v1));
124: } else {
125: PetscCall(VecRestoreArrayRead(Vec1, &v1));
126: PetscCall(VecRestoreArrayRead(Vec2, &v2));
127: }
128: }
129: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)Vec1), n_lt, lt, PETSC_OWN_POINTER, S));
130: PetscFunctionReturn(PETSC_SUCCESS);
131: }
133: /*@
134: VecWhichGreaterThan - Creates an index set containing the indices
135: where the vectors `Vec1` > `Vec2`
137: Collective
139: Input Parameters:
140: + Vec1 - the first vector to compare
141: - Vec2 - the second vector to compare
143: Output Parameter:
144: . S - The index set containing the indices i where vec1[i] > vec2[i]
146: Level: advanced
148: Notes:
149: The two vectors must have the same parallel layout
151: For complex numbers this only compares the real part
153: .seealso: `Vec`
154: @*/
155: PetscErrorCode VecWhichGreaterThan(Vec Vec1, Vec Vec2, IS *S)
156: {
157: PetscInt i, n_gt = 0;
158: PetscInt n, low, high;
159: PetscInt *gt = NULL;
160: const PetscScalar *v1, *v2;
162: PetscFunctionBegin;
165: PetscCheckSameComm(Vec1, 1, Vec2, 2);
166: VecCheckSameSize(Vec1, 1, Vec2, 2);
168: PetscCall(VecGetOwnershipRange(Vec1, &low, &high));
169: PetscCall(VecGetLocalSize(Vec1, &n));
170: if (n > 0) {
171: if (Vec1 == Vec2) {
172: PetscCall(VecGetArrayRead(Vec1, &v1));
173: v2 = v1;
174: } else {
175: PetscCall(VecGetArrayRead(Vec1, &v1));
176: PetscCall(VecGetArrayRead(Vec2, &v2));
177: }
179: PetscCall(PetscMalloc1(n, >));
181: for (i = 0; i < n; ++i) {
182: if (PetscRealPart(v1[i]) > PetscRealPart(v2[i])) {
183: gt[n_gt] = low + i;
184: ++n_gt;
185: }
186: }
188: if (Vec1 == Vec2) {
189: PetscCall(VecRestoreArrayRead(Vec1, &v1));
190: } else {
191: PetscCall(VecRestoreArrayRead(Vec1, &v1));
192: PetscCall(VecRestoreArrayRead(Vec2, &v2));
193: }
194: }
195: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)Vec1), n_gt, gt, PETSC_OWN_POINTER, S));
196: PetscFunctionReturn(PETSC_SUCCESS);
197: }
199: /*@
200: VecWhichBetween - Creates an index set containing the indices
201: where `VecLow` < `V` < `VecHigh`
203: Collective
205: Input Parameters:
206: + VecLow - lower bound
207: . V - Vector to compare
208: - VecHigh - higher bound
210: Output Parameter:
211: . S - The index set containing the indices i where veclow[i] < v[i] < vechigh[i]
213: Level: advanced
215: Notes:
216: The vectors must have the same parallel layout
218: For complex numbers this only compares the real part
220: .seealso: `Vec`
221: @*/
222: PetscErrorCode VecWhichBetween(Vec VecLow, Vec V, Vec VecHigh, IS *S)
223: {
224: PetscInt i, n_vm = 0;
225: PetscInt n, low, high;
226: PetscInt *vm = NULL;
227: const PetscScalar *v1, *v2, *vmiddle;
229: PetscFunctionBegin;
233: PetscCheckSameComm(V, 2, VecLow, 1);
234: PetscCheckSameComm(V, 2, VecHigh, 3);
235: VecCheckSameSize(V, 2, VecLow, 1);
236: VecCheckSameSize(V, 2, VecHigh, 3);
238: PetscCall(VecGetOwnershipRange(VecLow, &low, &high));
239: PetscCall(VecGetLocalSize(VecLow, &n));
240: if (n > 0) {
241: PetscCall(VecGetArrayRead(VecLow, &v1));
242: if (VecLow != VecHigh) {
243: PetscCall(VecGetArrayRead(VecHigh, &v2));
244: } else {
245: v2 = v1;
246: }
247: if (V != VecLow && V != VecHigh) {
248: PetscCall(VecGetArrayRead(V, &vmiddle));
249: } else if (V == VecLow) {
250: vmiddle = v1;
251: } else {
252: vmiddle = v2;
253: }
255: PetscCall(PetscMalloc1(n, &vm));
257: for (i = 0; i < n; ++i) {
258: if (PetscRealPart(v1[i]) < PetscRealPart(vmiddle[i]) && PetscRealPart(vmiddle[i]) < PetscRealPart(v2[i])) {
259: vm[n_vm] = low + i;
260: ++n_vm;
261: }
262: }
264: PetscCall(VecRestoreArrayRead(VecLow, &v1));
265: if (VecLow != VecHigh) PetscCall(VecRestoreArrayRead(VecHigh, &v2));
266: if (V != VecLow && V != VecHigh) PetscCall(VecRestoreArrayRead(V, &vmiddle));
267: }
268: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)V), n_vm, vm, PETSC_OWN_POINTER, S));
269: PetscFunctionReturn(PETSC_SUCCESS);
270: }
272: /*@
273: VecWhichBetweenOrEqual - Creates an index set containing the indices
274: where `VecLow` <= `V` <= `VecHigh`
276: Collective
278: Input Parameters:
279: + VecLow - lower bound
280: . V - Vector to compare
281: - VecHigh - higher bound
283: Output Parameter:
284: . S - The index set containing the indices i where veclow[i] <= v[i] <= vechigh[i]
286: Level: advanced
288: .seealso: `Vec`
289: @*/
290: PetscErrorCode VecWhichBetweenOrEqual(Vec VecLow, Vec V, Vec VecHigh, IS *S)
291: {
292: PetscInt i, n_vm = 0;
293: PetscInt n, low, high;
294: PetscInt *vm = NULL;
295: const PetscScalar *v1, *v2, *vmiddle;
297: PetscFunctionBegin;
301: PetscCheckSameComm(V, 2, VecLow, 1);
302: PetscCheckSameComm(V, 2, VecHigh, 3);
303: VecCheckSameSize(V, 2, VecLow, 1);
304: VecCheckSameSize(V, 2, VecHigh, 3);
306: PetscCall(VecGetOwnershipRange(VecLow, &low, &high));
307: PetscCall(VecGetLocalSize(VecLow, &n));
308: if (n > 0) {
309: PetscCall(VecGetArrayRead(VecLow, &v1));
310: if (VecLow != VecHigh) {
311: PetscCall(VecGetArrayRead(VecHigh, &v2));
312: } else {
313: v2 = v1;
314: }
315: if (V != VecLow && V != VecHigh) {
316: PetscCall(VecGetArrayRead(V, &vmiddle));
317: } else if (V == VecLow) {
318: vmiddle = v1;
319: } else {
320: vmiddle = v2;
321: }
323: PetscCall(PetscMalloc1(n, &vm));
325: for (i = 0; i < n; ++i) {
326: if (PetscRealPart(v1[i]) <= PetscRealPart(vmiddle[i]) && PetscRealPart(vmiddle[i]) <= PetscRealPart(v2[i])) {
327: vm[n_vm] = low + i;
328: ++n_vm;
329: }
330: }
332: PetscCall(VecRestoreArrayRead(VecLow, &v1));
333: if (VecLow != VecHigh) PetscCall(VecRestoreArrayRead(VecHigh, &v2));
334: if (V != VecLow && V != VecHigh) PetscCall(VecRestoreArrayRead(V, &vmiddle));
335: }
336: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)V), n_vm, vm, PETSC_OWN_POINTER, S));
337: PetscFunctionReturn(PETSC_SUCCESS);
338: }
340: /*@
341: VecWhichInactive - Creates an `IS` based on a set of vectors
343: Collective
345: Input Parameters:
346: + VecLow - lower bound
347: . V - Vector to compare
348: . D - Direction to compare
349: . VecHigh - higher bound
350: - Strong - indicator for applying strongly inactive test
352: Output Parameter:
353: . S - The index set containing the indices i where the bound is inactive
355: Level: advanced
357: Notes:
358: Creates an index set containing the indices where one of the following holds\:
359: .vb
360: - VecLow(i) < V(i) < VecHigh(i)
361: - VecLow(i) = V(i) and D(i) <= 0 (< 0 when Strong is true)
362: - VecHigh(i) = V(i) and D(i) >= 0 (> 0 when Strong is true)
363: .ve
365: .seealso: `Vec`
366: @*/
367: PetscErrorCode VecWhichInactive(Vec VecLow, Vec V, Vec D, Vec VecHigh, PetscBool Strong, IS *S)
368: {
369: PetscInt i, n_vm = 0;
370: PetscInt n, low, high;
371: PetscInt *vm = NULL;
372: const PetscScalar *v1, *v2, *v, *d;
374: PetscFunctionBegin;
375: if (!VecLow && !VecHigh) {
376: PetscCall(VecGetOwnershipRange(V, &low, &high));
377: PetscCall(ISCreateStride(PetscObjectComm((PetscObject)V), high - low, low, 1, S));
378: PetscFunctionReturn(PETSC_SUCCESS);
379: }
384: PetscCheckSameComm(V, 2, VecLow, 1);
385: PetscCheckSameComm(V, 2, D, 3);
386: PetscCheckSameComm(V, 2, VecHigh, 4);
387: VecCheckSameSize(V, 2, VecLow, 1);
388: VecCheckSameSize(V, 2, D, 3);
389: VecCheckSameSize(V, 2, VecHigh, 4);
391: PetscCall(VecGetOwnershipRange(VecLow, &low, &high));
392: PetscCall(VecGetLocalSize(VecLow, &n));
393: if (n > 0) {
394: PetscCall(VecGetArrayRead(VecLow, &v1));
395: if (VecLow != VecHigh) {
396: PetscCall(VecGetArrayRead(VecHigh, &v2));
397: } else {
398: v2 = v1;
399: }
400: if (V != VecLow && V != VecHigh) {
401: PetscCall(VecGetArrayRead(V, &v));
402: } else if (V == VecLow) {
403: v = v1;
404: } else {
405: v = v2;
406: }
407: if (D != VecLow && D != VecHigh && D != V) {
408: PetscCall(VecGetArrayRead(D, &d));
409: } else if (D == VecLow) {
410: d = v1;
411: } else if (D == VecHigh) {
412: d = v2;
413: } else {
414: d = v;
415: }
417: PetscCall(PetscMalloc1(n, &vm));
419: if (Strong) {
420: for (i = 0; i < n; ++i) {
421: if (PetscRealPart(v1[i]) < PetscRealPart(v[i]) && PetscRealPart(v[i]) < PetscRealPart(v2[i])) {
422: vm[n_vm] = low + i;
423: ++n_vm;
424: } else if (PetscRealPart(v1[i]) == PetscRealPart(v[i]) && PetscRealPart(d[i]) < 0) {
425: vm[n_vm] = low + i;
426: ++n_vm;
427: } else if (PetscRealPart(v2[i]) == PetscRealPart(v[i]) && PetscRealPart(d[i]) > 0) {
428: vm[n_vm] = low + i;
429: ++n_vm;
430: }
431: }
432: } else {
433: for (i = 0; i < n; ++i) {
434: if (PetscRealPart(v1[i]) < PetscRealPart(v[i]) && PetscRealPart(v[i]) < PetscRealPart(v2[i])) {
435: vm[n_vm] = low + i;
436: ++n_vm;
437: } else if (PetscRealPart(v1[i]) == PetscRealPart(v[i]) && PetscRealPart(d[i]) <= 0) {
438: vm[n_vm] = low + i;
439: ++n_vm;
440: } else if (PetscRealPart(v2[i]) == PetscRealPart(v[i]) && PetscRealPart(d[i]) >= 0) {
441: vm[n_vm] = low + i;
442: ++n_vm;
443: }
444: }
445: }
447: PetscCall(VecRestoreArrayRead(VecLow, &v1));
448: if (VecLow != VecHigh) PetscCall(VecRestoreArrayRead(VecHigh, &v2));
449: if (V != VecLow && V != VecHigh) PetscCall(VecRestoreArrayRead(V, &v));
450: if (D != VecLow && D != VecHigh && D != V) PetscCall(VecRestoreArrayRead(D, &d));
451: }
452: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)V), n_vm, vm, PETSC_OWN_POINTER, S));
453: PetscFunctionReturn(PETSC_SUCCESS);
454: }
456: /*@
457: VecISAXPY - Adds a reduced vector to the appropriate elements of a full-space vector.
458: vfull[is[i]] += alpha*vreduced[i]
460: Input Parameters:
461: + vfull - the full-space vector
462: . is - the index set for the reduced space
463: . alpha - the scalar coefficient
464: - vreduced - the reduced-space vector
466: Output Parameter:
467: . vfull - the sum of the full-space vector and reduced-space vector
469: Level: advanced
471: Notes:
472: The index set identifies entries in the global vector.
473: Negative indices are skipped; indices outside the ownership range of `vfull` will raise an error.
475: .seealso: `VecISCopy()`, `VecISSet()`, `VecAXPY()`
476: @*/
477: PetscErrorCode VecISAXPY(Vec vfull, IS is, PetscScalar alpha, Vec vreduced)
478: {
479: PetscInt nfull, nreduced;
480: PetscBool sorted = PETSC_FALSE;
482: PetscFunctionBegin;
486: PetscCall(VecGetSize(vfull, &nfull));
487: PetscCall(VecGetSize(vreduced, &nreduced));
488: if (nfull == nreduced) PetscCall(ISGetInfo(is, IS_SORTED, IS_GLOBAL, PETSC_TRUE, &sorted));
489: if (sorted) PetscCall(VecAXPY(vfull, alpha, vreduced));
490: else {
491: PetscScalar *y;
492: const PetscScalar *x;
493: PetscInt i, n, m, rstart, rend;
494: const PetscInt *id;
496: PetscCall(VecGetArray(vfull, &y));
497: PetscCall(VecGetArrayRead(vreduced, &x));
498: PetscCall(ISGetIndices(is, &id));
499: PetscCall(ISGetLocalSize(is, &n));
500: PetscCall(VecGetLocalSize(vreduced, &m));
501: PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_SUP, "IS local length %" PetscInt_FMT " not equal to Vec local length %" PetscInt_FMT, n, m);
502: PetscCall(VecGetOwnershipRange(vfull, &rstart, &rend));
503: y -= rstart;
504: if (alpha == 1.0) {
505: for (i = 0; i < n; ++i) {
506: if (id[i] < 0) continue;
507: PetscCheck(id[i] >= rstart && id[i] < rend, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only owned values supported");
508: y[id[i]] += x[i];
509: }
510: } else {
511: for (i = 0; i < n; ++i) {
512: if (id[i] < 0) continue;
513: PetscCheck(id[i] >= rstart && id[i] < rend, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only owned values supported");
514: y[id[i]] += alpha * x[i];
515: }
516: }
517: y += rstart;
518: PetscCall(ISRestoreIndices(is, &id));
519: PetscCall(VecRestoreArray(vfull, &y));
520: PetscCall(VecRestoreArrayRead(vreduced, &x));
521: }
522: PetscFunctionReturn(PETSC_SUCCESS);
523: }
525: /*@
526: VecISCopy - Copies between a reduced vector and the appropriate elements of a full-space vector.
528: Input Parameters:
529: + vfull - the full-space vector
530: . is - the index set for the reduced space
531: . mode - the direction of copying, `SCATTER_FORWARD` or `SCATTER_REVERSE`
532: - vreduced - the reduced-space vector
534: Output Parameter:
535: . vfull - the sum of the full-space vector and reduced-space vector
537: Level: advanced
539: Notes:
540: The index set identifies entries in the global vector.
541: Negative indices are skipped; indices outside the ownership range of `vfull` will raise an error.
542: .vb
543: mode == SCATTER_FORWARD: vfull[is[i]] = vreduced[i]
544: mode == SCATTER_REVERSE: vreduced[i] = vfull[is[i]]
545: .ve
547: .seealso: `VecISSet()`, `VecISAXPY()`, `VecCopy()`
548: @*/
549: PetscErrorCode VecISCopy(Vec vfull, IS is, ScatterMode mode, Vec vreduced)
550: {
551: PetscInt nfull, nreduced;
552: PetscBool sorted = PETSC_FALSE;
554: PetscFunctionBegin;
558: PetscCall(VecGetSize(vfull, &nfull));
559: PetscCall(VecGetSize(vreduced, &nreduced));
560: if (nfull == nreduced) PetscCall(ISGetInfo(is, IS_SORTED, IS_GLOBAL, PETSC_TRUE, &sorted));
561: if (sorted) {
562: if (mode == SCATTER_FORWARD) {
563: PetscCall(VecCopy(vreduced, vfull));
564: } else {
565: PetscCall(VecCopy(vfull, vreduced));
566: }
567: } else {
568: const PetscInt *id;
569: PetscInt i, n, m, rstart, rend;
571: PetscCall(ISGetIndices(is, &id));
572: PetscCall(ISGetLocalSize(is, &n));
573: PetscCall(VecGetLocalSize(vreduced, &m));
574: PetscCall(VecGetOwnershipRange(vfull, &rstart, &rend));
575: PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_SUP, "IS local length %" PetscInt_FMT " not equal to Vec local length %" PetscInt_FMT, n, m);
576: if (mode == SCATTER_FORWARD) {
577: PetscScalar *y;
578: const PetscScalar *x;
580: PetscCall(VecGetArray(vfull, &y));
581: PetscCall(VecGetArrayRead(vreduced, &x));
582: y -= rstart;
583: for (i = 0; i < n; ++i) {
584: if (id[i] < 0) continue;
585: PetscCheck(id[i] >= rstart && id[i] < rend, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only owned values supported");
586: y[id[i]] = x[i];
587: }
588: y += rstart;
589: PetscCall(VecRestoreArrayRead(vreduced, &x));
590: PetscCall(VecRestoreArray(vfull, &y));
591: } else if (mode == SCATTER_REVERSE) {
592: PetscScalar *x;
593: const PetscScalar *y;
595: PetscCall(VecGetArrayRead(vfull, &y));
596: PetscCall(VecGetArray(vreduced, &x));
597: for (i = 0; i < n; ++i) {
598: if (id[i] < 0) continue;
599: PetscCheck(id[i] >= rstart && id[i] < rend, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only owned values supported");
600: x[i] = y[id[i] - rstart];
601: }
602: PetscCall(VecRestoreArray(vreduced, &x));
603: PetscCall(VecRestoreArrayRead(vfull, &y));
604: } else SETERRQ(PetscObjectComm((PetscObject)vfull), PETSC_ERR_ARG_WRONG, "Only forward or reverse modes are legal");
605: PetscCall(ISRestoreIndices(is, &id));
606: }
607: PetscFunctionReturn(PETSC_SUCCESS);
608: }
610: /*@
611: ISComplementVec - Creates the complement of the index set relative to a layout defined by a `Vec`
613: Collective
615: Input Parameters:
616: + S - a PETSc `IS`
617: - V - the reference vector space
619: Output Parameter:
620: . T - the complement of S
622: Level: advanced
624: .seealso: `IS`, `Vec`, `ISCreateGeneral()`
625: @*/
626: PetscErrorCode ISComplementVec(IS S, Vec V, IS *T)
627: {
628: PetscInt start, end;
630: PetscFunctionBegin;
631: PetscCall(VecGetOwnershipRange(V, &start, &end));
632: PetscCall(ISComplement(S, start, end, T));
633: PetscFunctionReturn(PETSC_SUCCESS);
634: }
636: /*@
637: VecISSet - Sets the elements of a vector, specified by an index set, to a constant
639: Input Parameters:
640: + V - the vector
641: . S - index set for the locations in the vector
642: - c - the constant
644: Level: advanced
646: Notes:
647: The index set identifies entries in the global vector.
648: Negative indices are skipped; indices outside the ownership range of V will raise an error.
650: .seealso: `VecISCopy()`, `VecISAXPY()`, `VecSet()`
651: @*/
652: PetscErrorCode VecISSet(Vec V, IS S, PetscScalar c)
653: {
654: PetscInt nloc, low, high, i;
655: const PetscInt *s;
656: PetscScalar *v;
658: PetscFunctionBegin;
659: if (!S) PetscFunctionReturn(PETSC_SUCCESS); /* simply return with no-op if the index set is NULL */
664: PetscCall(VecGetOwnershipRange(V, &low, &high));
665: PetscCall(ISGetLocalSize(S, &nloc));
666: PetscCall(ISGetIndices(S, &s));
667: PetscCall(VecGetArray(V, &v));
668: for (i = 0; i < nloc; ++i) {
669: if (s[i] < 0) continue;
670: PetscCheck(s[i] >= low && s[i] < high, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only owned values supported");
671: v[s[i] - low] = c;
672: }
673: PetscCall(ISRestoreIndices(S, &s));
674: PetscCall(VecRestoreArray(V, &v));
675: PetscFunctionReturn(PETSC_SUCCESS);
676: }
678: #if !defined(PETSC_USE_COMPLEX)
679: /*@C
680: VecBoundGradientProjection - Projects vector according to this definition.
681: If XL[i] < X[i] < XU[i], then GP[i] = G[i];
682: If X[i] <= XL[i], then GP[i] = min(G[i],0);
683: If X[i] >= XU[i], then GP[i] = max(G[i],0);
685: Input Parameters:
686: + G - current gradient vector
687: . X - current solution vector with XL[i] <= X[i] <= XU[i]
688: . XL - lower bounds
689: - XU - upper bounds
691: Output Parameter:
692: . GP - gradient projection vector
694: Level: advanced
696: Note:
697: `GP` may be the same vector as `G`
699: .seealso: `Vec`
700: @*/
701: PetscErrorCode VecBoundGradientProjection(Vec G, Vec X, Vec XL, Vec XU, Vec GP)
702: {
703: PetscInt n, i;
704: const PetscReal *xptr, *xlptr, *xuptr;
705: PetscReal *gptr, *gpptr;
706: PetscReal xval, gpval;
708: /* Project variables at the lower and upper bound */
709: PetscFunctionBegin;
715: if (!XL && !XU) {
716: PetscCall(VecCopy(G, GP));
717: PetscFunctionReturn(PETSC_SUCCESS);
718: }
720: PetscCall(VecGetLocalSize(X, &n));
722: PetscCall(VecGetArrayRead(X, &xptr));
723: PetscCall(VecGetArrayRead(XL, &xlptr));
724: PetscCall(VecGetArrayRead(XU, &xuptr));
725: PetscCall(VecGetArrayPair(G, GP, &gptr, &gpptr));
727: for (i = 0; i < n; ++i) {
728: gpval = gptr[i];
729: xval = xptr[i];
730: if (gpval > 0.0 && xval <= xlptr[i]) {
731: gpval = 0.0;
732: } else if (gpval < 0.0 && xval >= xuptr[i]) {
733: gpval = 0.0;
734: }
735: gpptr[i] = gpval;
736: }
738: PetscCall(VecRestoreArrayRead(X, &xptr));
739: PetscCall(VecRestoreArrayRead(XL, &xlptr));
740: PetscCall(VecRestoreArrayRead(XU, &xuptr));
741: PetscCall(VecRestoreArrayPair(G, GP, &gptr, &gpptr));
742: PetscFunctionReturn(PETSC_SUCCESS);
743: }
744: #endif
746: /*@
747: VecStepMaxBounded - See below
749: Collective
751: Input Parameters:
752: + X - vector with no negative entries
753: . XL - lower bounds
754: . XU - upper bounds
755: - DX - step direction, can have negative, positive or zero entries
757: Output Parameter:
758: . stepmax - minimum value so that X[i] + stepmax*DX[i] <= XL[i] or XU[i] <= X[i] + stepmax*DX[i]
760: Level: intermediate
762: .seealso: `Vec`
763: @*/
764: PetscErrorCode VecStepMaxBounded(Vec X, Vec DX, Vec XL, Vec XU, PetscReal *stepmax)
765: {
766: PetscInt i, nn;
767: const PetscScalar *xx, *dx, *xl, *xu;
768: PetscReal localmax = 0;
770: PetscFunctionBegin;
776: PetscCall(VecGetArrayRead(X, &xx));
777: PetscCall(VecGetArrayRead(XL, &xl));
778: PetscCall(VecGetArrayRead(XU, &xu));
779: PetscCall(VecGetArrayRead(DX, &dx));
780: PetscCall(VecGetLocalSize(X, &nn));
781: for (i = 0; i < nn; i++) {
782: if (PetscRealPart(dx[i]) > 0) {
783: localmax = PetscMax(localmax, PetscRealPart((xu[i] - xx[i]) / dx[i]));
784: } else if (PetscRealPart(dx[i]) < 0) {
785: localmax = PetscMax(localmax, PetscRealPart((xl[i] - xx[i]) / dx[i]));
786: }
787: }
788: PetscCall(VecRestoreArrayRead(X, &xx));
789: PetscCall(VecRestoreArrayRead(XL, &xl));
790: PetscCall(VecRestoreArrayRead(XU, &xu));
791: PetscCall(VecRestoreArrayRead(DX, &dx));
792: PetscCall(MPIU_Allreduce(&localmax, stepmax, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)X)));
793: PetscFunctionReturn(PETSC_SUCCESS);
794: }
796: /*@
797: VecStepBoundInfo - See below
799: Collective
801: Input Parameters:
802: + X - vector with no negative entries
803: . XL - lower bounds
804: . XU - upper bounds
805: - DX - step direction, can have negative, positive or zero entries
807: Output Parameters:
808: + boundmin - (may be `NULL` this it is not computed) maximum value so that XL[i] <= X[i] + boundmax*DX[i] <= XU[i]
809: . wolfemin - (may be `NULL` this it is not computed)
810: - 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]
812: Level: advanced
814: Note:
815: For complex numbers only compares the real part
817: .seealso: `Vec`
818: @*/
819: PetscErrorCode VecStepBoundInfo(Vec X, Vec DX, Vec XL, Vec XU, PetscReal *boundmin, PetscReal *wolfemin, PetscReal *boundmax)
820: {
821: PetscInt n, i;
822: const PetscScalar *x, *xl, *xu, *dx;
823: PetscReal t;
824: PetscReal localmin = PETSC_INFINITY, localwolfemin = PETSC_INFINITY, localmax = -1;
825: MPI_Comm comm;
827: PetscFunctionBegin;
833: PetscCall(VecGetArrayRead(X, &x));
834: PetscCall(VecGetArrayRead(XL, &xl));
835: PetscCall(VecGetArrayRead(XU, &xu));
836: PetscCall(VecGetArrayRead(DX, &dx));
837: PetscCall(VecGetLocalSize(X, &n));
838: for (i = 0; i < n; ++i) {
839: if (PetscRealPart(dx[i]) > 0 && PetscRealPart(xu[i]) < PETSC_INFINITY) {
840: t = PetscRealPart((xu[i] - x[i]) / dx[i]);
841: localmin = PetscMin(t, localmin);
842: if (localmin > 0) localwolfemin = PetscMin(t, localwolfemin);
843: localmax = PetscMax(t, localmax);
844: } else if (PetscRealPart(dx[i]) < 0 && PetscRealPart(xl[i]) > PETSC_NINFINITY) {
845: t = PetscRealPart((xl[i] - x[i]) / dx[i]);
846: localmin = PetscMin(t, localmin);
847: if (localmin > 0) localwolfemin = PetscMin(t, localwolfemin);
848: localmax = PetscMax(t, localmax);
849: }
850: }
852: PetscCall(VecRestoreArrayRead(X, &x));
853: PetscCall(VecRestoreArrayRead(XL, &xl));
854: PetscCall(VecRestoreArrayRead(XU, &xu));
855: PetscCall(VecRestoreArrayRead(DX, &dx));
856: PetscCall(PetscObjectGetComm((PetscObject)X, &comm));
858: if (boundmin) {
859: PetscCall(MPIU_Allreduce(&localmin, boundmin, 1, MPIU_REAL, MPIU_MIN, comm));
860: PetscCall(PetscInfo(X, "Step Bound Info: Closest Bound: %20.19e\n", (double)*boundmin));
861: }
862: if (wolfemin) {
863: PetscCall(MPIU_Allreduce(&localwolfemin, wolfemin, 1, MPIU_REAL, MPIU_MIN, comm));
864: PetscCall(PetscInfo(X, "Step Bound Info: Wolfe: %20.19e\n", (double)*wolfemin));
865: }
866: if (boundmax) {
867: PetscCall(MPIU_Allreduce(&localmax, boundmax, 1, MPIU_REAL, MPIU_MAX, comm));
868: if (*boundmax < 0) *boundmax = PETSC_INFINITY;
869: PetscCall(PetscInfo(X, "Step Bound Info: Max: %20.19e\n", (double)*boundmax));
870: }
871: PetscFunctionReturn(PETSC_SUCCESS);
872: }
874: /*@
875: VecStepMax - Returns the largest value so that x[i] + step*DX[i] >= 0 for all i
877: Collective
879: Input Parameters:
880: + X - vector with no negative entries
881: - DX - a step direction, can have negative, positive or zero entries
883: Output Parameter:
884: . step - largest value such that x[i] + step*DX[i] >= 0 for all i
886: Level: advanced
888: Note:
889: For complex numbers only compares the real part
891: .seealso: `Vec`
892: @*/
893: PetscErrorCode VecStepMax(Vec X, Vec DX, PetscReal *step)
894: {
895: PetscInt i, nn;
896: PetscReal stepmax = PETSC_INFINITY;
897: const PetscScalar *xx, *dx;
899: PetscFunctionBegin;
903: PetscCall(VecGetLocalSize(X, &nn));
904: PetscCall(VecGetArrayRead(X, &xx));
905: PetscCall(VecGetArrayRead(DX, &dx));
906: for (i = 0; i < nn; ++i) {
907: PetscCheck(PetscRealPart(xx[i]) >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vector must be positive");
908: if (PetscRealPart(dx[i]) < 0) stepmax = PetscMin(stepmax, PetscRealPart(-xx[i] / dx[i]));
909: }
910: PetscCall(VecRestoreArrayRead(X, &xx));
911: PetscCall(VecRestoreArrayRead(DX, &dx));
912: PetscCall(MPIU_Allreduce(&stepmax, step, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)X)));
913: PetscFunctionReturn(PETSC_SUCCESS);
914: }
916: /*@
917: VecPow - Replaces each component of a vector by x_i^p
919: Logically Collective
921: Input Parameters:
922: + v - the vector
923: - p - the exponent to use on each element
925: Level: intermediate
927: .seealso: `Vec`
928: @*/
929: PetscErrorCode VecPow(Vec v, PetscScalar p)
930: {
931: PetscInt n, i;
932: PetscScalar *v1;
934: PetscFunctionBegin;
937: PetscCall(VecGetArray(v, &v1));
938: PetscCall(VecGetLocalSize(v, &n));
940: if (1.0 == p) {
941: } else if (-1.0 == p) {
942: for (i = 0; i < n; ++i) v1[i] = 1.0 / v1[i];
943: } else if (0.0 == p) {
944: for (i = 0; i < n; ++i) {
945: /* Not-a-number left alone
946: Infinity set to one */
947: if (v1[i] == v1[i]) v1[i] = 1.0;
948: }
949: } else if (0.5 == p) {
950: for (i = 0; i < n; ++i) {
951: if (PetscRealPart(v1[i]) >= 0) {
952: v1[i] = PetscSqrtScalar(v1[i]);
953: } else {
954: v1[i] = PETSC_INFINITY;
955: }
956: }
957: } else if (-0.5 == p) {
958: for (i = 0; i < n; ++i) {
959: if (PetscRealPart(v1[i]) >= 0) {
960: v1[i] = 1.0 / PetscSqrtScalar(v1[i]);
961: } else {
962: v1[i] = PETSC_INFINITY;
963: }
964: }
965: } else if (2.0 == p) {
966: for (i = 0; i < n; ++i) v1[i] *= v1[i];
967: } else if (-2.0 == p) {
968: for (i = 0; i < n; ++i) v1[i] = 1.0 / (v1[i] * v1[i]);
969: } else {
970: for (i = 0; i < n; ++i) {
971: if (PetscRealPart(v1[i]) >= 0) {
972: v1[i] = PetscPowScalar(v1[i], p);
973: } else {
974: v1[i] = PETSC_INFINITY;
975: }
976: }
977: }
978: PetscCall(VecRestoreArray(v, &v1));
979: PetscFunctionReturn(PETSC_SUCCESS);
980: }
982: /*@
983: VecMedian - Computes the componentwise median of three vectors
984: and stores the result in this vector. Used primarily for projecting
985: a vector within upper and lower bounds.
987: Logically Collective
989: Input Parameters:
990: + Vec1 - The first vector
991: . Vec2 - The second vector
992: - Vec3 - The third vector
994: Output Parameter:
995: . VMedian - The median vector (this can be any one of the input vectors)
997: Level: advanced
999: .seealso: `Vec`
1000: @*/
1001: PetscErrorCode VecMedian(Vec Vec1, Vec Vec2, Vec Vec3, Vec VMedian)
1002: {
1003: PetscInt i, n, low1, high1;
1004: const PetscScalar *v1, *v2, *v3;
1005: PetscScalar *vmed;
1007: PetscFunctionBegin;
1013: if (!Vec1 && !Vec3) {
1014: PetscCall(VecCopy(Vec2, VMedian));
1015: PetscFunctionReturn(PETSC_SUCCESS);
1016: }
1017: if (Vec1 == Vec2 || Vec1 == Vec3) {
1018: PetscCall(VecCopy(Vec1, VMedian));
1019: PetscFunctionReturn(PETSC_SUCCESS);
1020: }
1021: if (Vec2 == Vec3) {
1022: PetscCall(VecCopy(Vec2, VMedian));
1023: PetscFunctionReturn(PETSC_SUCCESS);
1024: }
1026: /* Assert that Vec1 != Vec2 and Vec2 != Vec3 */
1031: PetscCheckSameType(Vec1, 1, Vec2, 2);
1032: PetscCheckSameType(Vec1, 1, Vec3, 3);
1033: PetscCheckSameType(Vec1, 1, VMedian, 4);
1034: PetscCheckSameComm(Vec1, 1, Vec2, 2);
1035: PetscCheckSameComm(Vec1, 1, Vec3, 3);
1036: PetscCheckSameComm(Vec1, 1, VMedian, 4);
1037: VecCheckSameSize(Vec1, 1, Vec2, 2);
1038: VecCheckSameSize(Vec1, 1, Vec3, 3);
1039: VecCheckSameSize(Vec1, 1, VMedian, 4);
1041: PetscCall(VecGetOwnershipRange(Vec1, &low1, &high1));
1042: PetscCall(VecGetLocalSize(Vec1, &n));
1043: if (n > 0) {
1044: PetscCall(VecGetArray(VMedian, &vmed));
1045: if (Vec1 != VMedian) {
1046: PetscCall(VecGetArrayRead(Vec1, &v1));
1047: } else {
1048: v1 = vmed;
1049: }
1050: if (Vec2 != VMedian) {
1051: PetscCall(VecGetArrayRead(Vec2, &v2));
1052: } else {
1053: v2 = vmed;
1054: }
1055: if (Vec3 != VMedian) {
1056: PetscCall(VecGetArrayRead(Vec3, &v3));
1057: } else {
1058: v3 = vmed;
1059: }
1061: 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])));
1063: PetscCall(VecRestoreArray(VMedian, &vmed));
1064: if (VMedian != Vec1) PetscCall(VecRestoreArrayRead(Vec1, &v1));
1065: if (VMedian != Vec2) PetscCall(VecRestoreArrayRead(Vec2, &v2));
1066: if (VMedian != Vec3) PetscCall(VecRestoreArrayRead(Vec3, &v3));
1067: }
1068: PetscFunctionReturn(PETSC_SUCCESS);
1069: }