Actual source code: relative.c


  2: #include <petsc/private/vecimpl.h>
  3: #include "../src/vec/vec/utils/tagger/impls/simple.h"

  5: static PetscErrorCode VecTaggerComputeBoxes_Relative(VecTagger tagger, Vec vec, PetscInt *numBoxes, VecTaggerBox **boxes, PetscBool *listed)
  6: {
  7:   VecTagger_Simple  *smpl = (VecTagger_Simple *)tagger->data;
  8:   PetscInt           bs, i, j, k, n;
  9:   VecTaggerBox      *bxs;
 10:   const PetscScalar *vArray;

 12:   VecTaggerGetBlockSize(tagger, &bs);
 13:   *numBoxes = 1;
 14:   PetscMalloc1(bs, &bxs);
 15:   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:   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:   VecRestoreArrayRead(vec, &vArray);
 40:   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:   return 0;
 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:   VecTaggerSetBox_Simple(tagger, box);
 73:   return 0;
 74: }

 76: /*@C
 77:   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.

 79:   Logically Collective

 81:   Input Parameter:
 82: . tagger - the VecTagger context

 84:   Output Parameter:
 85: . box - a blocksize list of VecTaggerBox boxes

 87:   Level: advanced

 89: .seealso: `VecTaggerRelativeSetBox()`
 90: @*/
 91: PetscErrorCode VecTaggerRelativeGetBox(VecTagger tagger, const VecTaggerBox **box)
 92: {
 93:   VecTaggerGetBox_Simple(tagger, box);
 94:   return 0;
 95: }

 97: PETSC_INTERN PetscErrorCode VecTaggerCreate_Relative(VecTagger tagger)
 98: {
 99:   VecTaggerCreate_Simple(tagger);
100:   tagger->ops->computeboxes = VecTaggerComputeBoxes_Relative;
101:   return 0;
102: }