15c6c1daeSBarry Smith 2999739cfSJacob Faibussowitsch #include <petsc/private/drawimpl.h> /*I "petscdraw.h" I*/ 35c6c1daeSBarry Smith 45c6c1daeSBarry Smith /* 55c6c1daeSBarry Smith val is the label value. sep is the separation to the next (or previous) 65c6c1daeSBarry Smith label; this is useful in determining how many significant figures to 75c6c1daeSBarry Smith keep. 85c6c1daeSBarry Smith */ 9d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscADefLabel(PetscReal val, PetscReal sep, char **p) 10d71ae5a4SJacob Faibussowitsch { 115c6c1daeSBarry Smith static char buf[40]; 125c6c1daeSBarry Smith 135c6c1daeSBarry Smith PetscFunctionBegin; 145c6c1daeSBarry Smith /* Find the string */ 15ba129bcfSBarry Smith if (PetscAbsReal(val) / sep < 1.e-4) { 169371c9d4SSatish Balay buf[0] = '0'; 179371c9d4SSatish Balay buf[1] = 0; 1836c9bc0dSBarry Smith } else { 19a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, PETSC_STATIC_ARRAY_LENGTH(buf), "%0.1e", (double)val)); 209566063dSJacob Faibussowitsch PetscCall(PetscStripZerosPlus(buf)); 219566063dSJacob Faibussowitsch PetscCall(PetscStripe0(buf)); 229566063dSJacob Faibussowitsch PetscCall(PetscStripInitialZero(buf)); 239566063dSJacob Faibussowitsch PetscCall(PetscStripAllZeros(buf)); 249566063dSJacob Faibussowitsch PetscCall(PetscStripTrailingZeros(buf)); 255c6c1daeSBarry Smith } 265c6c1daeSBarry Smith *p = buf; 27*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 285c6c1daeSBarry Smith } 295c6c1daeSBarry Smith 305c6c1daeSBarry Smith /* Finds "nice" locations for the ticks */ 31d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscADefTicks(PetscReal low, PetscReal high, int num, int *ntick, PetscReal *tickloc, int maxtick) 32d71ae5a4SJacob Faibussowitsch { 335c6c1daeSBarry Smith int i, power; 34e5ab1681SLisandro Dalcin PetscReal x = 0.0, base = 0.0, eps; 355c6c1daeSBarry Smith 365c6c1daeSBarry Smith PetscFunctionBegin; 379566063dSJacob Faibussowitsch PetscCall(PetscAGetBase(low, high, num, &base, &power)); 389566063dSJacob Faibussowitsch PetscCall(PetscAGetNice(low, base, -1, &x)); 395c6c1daeSBarry Smith 405c6c1daeSBarry Smith /* Values are of the form j * base */ 415c6c1daeSBarry Smith /* Find the starting value */ 425c6c1daeSBarry Smith if (x < low) x += base; 435c6c1daeSBarry Smith 449371c9d4SSatish Balay i = 0; 459371c9d4SSatish Balay eps = base / 10; 46e5ab1681SLisandro Dalcin while (i < maxtick && x <= high + eps) { 479371c9d4SSatish Balay tickloc[i++] = x; 489371c9d4SSatish Balay x += base; 495c6c1daeSBarry Smith } 505c6c1daeSBarry Smith *ntick = i; 51e5ab1681SLisandro Dalcin tickloc[i - 1] = PetscMin(tickloc[i - 1], high); 525c6c1daeSBarry Smith 5348a46eb9SPierre Jolivet if (i < 2 && num < 10) PetscCall(PetscADefTicks(low, high, num + 1, ntick, tickloc, maxtick)); 54*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 555c6c1daeSBarry Smith } 565c6c1daeSBarry Smith 575c6c1daeSBarry Smith #define EPS 1.e-6 585c6c1daeSBarry Smith 59d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscExp10(PetscReal d, PetscReal *result) 60d71ae5a4SJacob Faibussowitsch { 615c6c1daeSBarry Smith PetscFunctionBegin; 6285ec1a3cSBarry Smith *result = PetscPowReal((PetscReal)10.0, d); 63*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 645c6c1daeSBarry Smith } 655c6c1daeSBarry Smith 66d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscMod(PetscReal x, PetscReal y, PetscReal *result) 67d71ae5a4SJacob Faibussowitsch { 685c6c1daeSBarry Smith int i; 695c6c1daeSBarry Smith 705c6c1daeSBarry Smith PetscFunctionBegin; 715c6c1daeSBarry Smith if (y == 1) { 725c6c1daeSBarry Smith *result = 0.0; 73*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 745c6c1daeSBarry Smith } 755c6c1daeSBarry Smith i = ((int)x) / ((int)y); 765c6c1daeSBarry Smith x = x - i * y; 775c6c1daeSBarry Smith while (x > y) x -= y; 785c6c1daeSBarry Smith *result = x; 79*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 805c6c1daeSBarry Smith } 815c6c1daeSBarry Smith 82d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCopysign(PetscReal a, PetscReal b, PetscReal *result) 83d71ae5a4SJacob Faibussowitsch { 845c6c1daeSBarry Smith PetscFunctionBegin; 855c6c1daeSBarry Smith if (b >= 0) *result = a; 865c6c1daeSBarry Smith else *result = -a; 87*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 885c6c1daeSBarry Smith } 895c6c1daeSBarry Smith 905c6c1daeSBarry Smith /* 915c6c1daeSBarry Smith Given a value "in" and a "base", return a nice value. 925c6c1daeSBarry Smith based on "sign", extend up (+1) or down (-1) 935c6c1daeSBarry Smith */ 94d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAGetNice(PetscReal in, PetscReal base, int sign, PetscReal *result) 95d71ae5a4SJacob Faibussowitsch { 965c6c1daeSBarry Smith PetscReal etmp, s, s2, m; 975c6c1daeSBarry Smith 985c6c1daeSBarry Smith PetscFunctionBegin; 999566063dSJacob Faibussowitsch PetscCall(PetscCopysign(0.5, (double)sign, &s)); 1005c6c1daeSBarry Smith etmp = in / base + 0.5 + s; 1019566063dSJacob Faibussowitsch PetscCall(PetscCopysign(0.5, etmp, &s)); 1029566063dSJacob Faibussowitsch PetscCall(PetscCopysign(EPS * etmp, (double)sign, &s2)); 1035c6c1daeSBarry Smith etmp = etmp - 0.5 + s - s2; 1049566063dSJacob Faibussowitsch PetscCall(PetscMod(etmp, 1.0, &m)); 1055c6c1daeSBarry Smith etmp = base * (etmp - m); 1065c6c1daeSBarry Smith *result = etmp; 107*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1085c6c1daeSBarry Smith } 1095c6c1daeSBarry Smith 110d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAGetBase(PetscReal vmin, PetscReal vmax, int num, PetscReal *Base, int *power) 111d71ae5a4SJacob Faibussowitsch { 1125c6c1daeSBarry Smith PetscReal base, ftemp, e10; 1135c6c1daeSBarry Smith static PetscReal base_try[5] = {10.0, 5.0, 2.0, 1.0, 0.5}; 1145c6c1daeSBarry Smith int i; 1155c6c1daeSBarry Smith 1165c6c1daeSBarry Smith PetscFunctionBegin; 1175c6c1daeSBarry Smith /* labels of the form n * BASE */ 1185c6c1daeSBarry Smith /* get an approximate value for BASE */ 1195c6c1daeSBarry Smith base = (vmax - vmin) / (double)(num + 1); 1205c6c1daeSBarry Smith 1215c6c1daeSBarry Smith /* make it of form m x 10^power, m in [1.0, 10) */ 1225c6c1daeSBarry Smith if (base <= 0.0) { 1235c6c1daeSBarry Smith base = PetscAbsReal(vmin); 1245c6c1daeSBarry Smith if (base < 1.0) base = 1.0; 1255c6c1daeSBarry Smith } 12677b4d14cSPeter Brune ftemp = PetscLog10Real((1.0 + EPS) * base); 1275c6c1daeSBarry Smith if (ftemp < 0.0) ftemp -= 1.0; 1285c6c1daeSBarry Smith *power = (int)ftemp; 1299566063dSJacob Faibussowitsch PetscCall(PetscExp10((double)-*power, &e10)); 1305c6c1daeSBarry Smith base = base * e10; 1315c6c1daeSBarry Smith if (base < 1.0) base = 1.0; 1325c6c1daeSBarry Smith /* now reduce it to one of 1, 2, or 5 */ 1335c6c1daeSBarry Smith for (i = 1; i < 5; i++) { 1345c6c1daeSBarry Smith if (base >= base_try[i]) { 1359566063dSJacob Faibussowitsch PetscCall(PetscExp10((double)*power, &e10)); 1365c6c1daeSBarry Smith base = base_try[i - 1] * e10; 1375c6c1daeSBarry Smith if (i == 1) *power = *power + 1; 1385c6c1daeSBarry Smith break; 1395c6c1daeSBarry Smith } 1405c6c1daeSBarry Smith } 1415c6c1daeSBarry Smith *Base = base; 142*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1435c6c1daeSBarry Smith } 144