Actual source code: relative.c
1: #include <petsc/private/vecimpl.h>
2: #include "../src/vec/vec/utils/tagger/impls/simple.h"
4: static PetscErrorCode VecTaggerComputeBoxes_Relative(VecTagger tagger, Vec vec, PetscInt *numBoxes, VecTaggerBox **boxes, PetscBool *listed)
5: {
6: VecTagger_Simple *smpl = (VecTagger_Simple *)tagger->data;
7: PetscInt bs, i, j, k, n;
8: VecTaggerBox *bxs;
9: const PetscScalar *vArray;
11: PetscFunctionBegin;
12: PetscCall(VecTaggerGetBlockSize(tagger, &bs));
13: *numBoxes = 1;
14: PetscCall(PetscMalloc1(bs, &bxs));
15: PetscCall(VecGetLocalSize(vec, &n));
16: n /= bs;
17: for (i = 0; i < bs; i++) {
18: #if !defined(PETSC_USE_COMPLEX)
19: bxs[i].min = PETSC_MAX_REAL;
20: bxs[i].max = PETSC_MIN_REAL;
21: #else
22: bxs[i].min = PetscCMPLX(PETSC_MAX_REAL, PETSC_MAX_REAL);
23: bxs[i].max = PetscCMPLX(PETSC_MIN_REAL, PETSC_MIN_REAL);
24: #endif
25: }
26: PetscCall(VecGetArrayRead(vec, &vArray));
27: for (i = 0, k = 0; i < n; i++) {
28: for (j = 0; j < bs; j++, k++) {
29: #if !defined(PETSC_USE_COMPLEX)
30: bxs[j].min = PetscMin(bxs[j].min, vArray[k]);
31: bxs[j].max = PetscMax(bxs[j].max, vArray[k]);
32: #else
33: bxs[j].min = PetscCMPLX(PetscMin(PetscRealPart(bxs[j].min), PetscRealPart(vArray[k])), PetscMin(PetscImaginaryPart(bxs[j].min), PetscImaginaryPart(vArray[k])));
34: bxs[j].max = PetscCMPLX(PetscMax(PetscRealPart(bxs[j].max), PetscRealPart(vArray[k])), PetscMax(PetscImaginaryPart(bxs[j].max), PetscImaginaryPart(vArray[k])));
35: #endif
36: }
37: }
38: for (i = 0; i < bs; i++) bxs[i].max = -bxs[i].max;
39: PetscCall(VecRestoreArrayRead(vec, &vArray));
40: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, (PetscReal *)bxs, 2 * (sizeof(PetscScalar) / sizeof(PetscReal)) * bs, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)tagger)));
41: for (i = 0; i < bs; i++) {
42: PetscScalar mins = bxs[i].min;
43: PetscScalar difs = -bxs[i].max - mins;
44: #if !defined(PETSC_USE_COMPLEX)
45: bxs[i].min = mins + smpl->box[i].min * difs;
46: bxs[i].max = mins + smpl->box[i].max * difs;
47: #else
48: bxs[i].min = mins + PetscCMPLX(PetscRealPart(smpl->box[i].min) * PetscRealPart(difs), PetscImaginaryPart(smpl->box[i].min) * PetscImaginaryPart(difs));
49: bxs[i].max = mins + PetscCMPLX(PetscRealPart(smpl->box[i].max) * PetscRealPart(difs), PetscImaginaryPart(smpl->box[i].max) * PetscImaginaryPart(difs));
50: #endif
51: }
52: *boxes = bxs;
53: if (listed) *listed = PETSC_TRUE;
54: PetscFunctionReturn(PETSC_SUCCESS);
55: }
57: /*@C
58: VecTaggerRelativeSetBox - Set the relative box defining the values to be tagged by the tagger, where relative boxes are subsets of [0,1] (or [0,1]+[0,1]i for complex scalars), where 0 indicates the smallest value present in the vector and 1 indicates the largest.
60: Logically Collective
62: Input Parameters:
63: + tagger - the VecTagger context
64: - box - a blocksize list of VecTaggerBox boxes
66: Level: advanced
68: .seealso: `VecTaggerRelativeGetBox()`
69: @*/
70: PetscErrorCode VecTaggerRelativeSetBox(VecTagger tagger, VecTaggerBox *box)
71: {
72: PetscFunctionBegin;
73: PetscCall(VecTaggerSetBox_Simple(tagger, box));
74: PetscFunctionReturn(PETSC_SUCCESS);
75: }
77: /*@C
78: VecTaggerRelativeGetBox - Get the relative box defining the values to be tagged by the tagger, where relative boxess are subsets of [0,1] (or [0,1]+[0,1]i for complex scalars), where 0 indicates the smallest value present in the vector and 1 indicates the largest.
80: Logically Collective
82: Input Parameter:
83: . tagger - the VecTagger context
85: Output Parameter:
86: . box - a blocksize list of VecTaggerBox boxes
88: Level: advanced
90: .seealso: `VecTaggerRelativeSetBox()`
91: @*/
92: PetscErrorCode VecTaggerRelativeGetBox(VecTagger tagger, const VecTaggerBox **box)
93: {
94: PetscFunctionBegin;
95: PetscCall(VecTaggerGetBox_Simple(tagger, box));
96: PetscFunctionReturn(PETSC_SUCCESS);
97: }
99: PETSC_INTERN PetscErrorCode VecTaggerCreate_Relative(VecTagger tagger)
100: {
101: PetscFunctionBegin;
102: PetscCall(VecTaggerCreate_Simple(tagger));
103: tagger->ops->computeboxes = VecTaggerComputeBoxes_Relative;
104: PetscFunctionReturn(PETSC_SUCCESS);
105: }