triangle_group_limit_set/draw.c

1595 lines
45 KiB
C

#include "main.h"
#define FMOD(x,y) (fmod(x,y) < 0 ? fmod(x,y) + y : fmod(x,y))
#define ANGLE_DIFF(x,y) (FMOD((x)-(y), 2*M_PI))
#define ANGLE_IN_INTERVAL(a,b,x) (ANGLE_DIFF(x,a) < ANGLE_DIFF(b,a))
#define FLIP(x,y) do {double tmp = x; x = y; y = tmp;} while(0)
// level 0: helper functions
int isInsideBB(DrawingContext *ctx, point_t p)
{
cairo_user_to_device(ctx->cairo, &p.x, &p.y);
return -p.x < ctx->dim->width && p.x < 3*ctx->dim->width && -p.y < ctx->dim->height && p.y < 3*ctx->dim->height;
}
vector_t cross(vector_t a, vector_t b)
{
vector_t result;
result.x[0] = a.x[1]*b.x[2] - a.x[2]*b.x[1];
result.x[1] = a.x[2]*b.x[0] - a.x[0]*b.x[2];
result.x[2] = a.x[0]*b.x[1] - a.x[1]*b.x[0];
return result;
}
vector_t apply(gsl_matrix *m, vector_t x)
{
vector_t out;
LOOP(i) out.x[i] = 0.0;
LOOP(i) LOOP(j) out.x[i] += gsl_matrix_get(m, i, j) * x.x[j];
return out;
}
vector_t apply_transpose(gsl_matrix *m, vector_t x)
{
vector_t out;
LOOP(i) out.x[i] = 0.0;
LOOP(i) LOOP(j) out.x[i] += gsl_matrix_get(m, j, i) * x.x[j];
return out;
}
int fixedPoints(DrawingContext *ctx, const char *word, vector_t *out)
{
gsl_matrix *tmp = getTempMatrix(ctx->ws);
gsl_matrix *ev = getTempMatrix(ctx->ws);
gsl_matrix **gen = getTempMatrices(ctx->ws, 6);
initializeTriangleGeneratorsCurrent(gen, ctx);
gsl_matrix_set_identity(tmp);
for(int i = 0; i < strlen(word); i++) {
if(word[i] >= 'a' && word[i] <= 'c')
multiply_right(tmp, gen[word[i]-'a'], ctx->ws);
else if(word[i] >= 'A' && word[i] <= 'C')
multiply_right(tmp, gen[word[i]-'A'+3], ctx->ws);
}
int count = real_eigenvectors(tmp, ev, ctx->ws);
LOOP(i) LOOP(j) out[i].x[j] = gsl_matrix_get(ev, j, i);
releaseTempMatrices(ctx->ws, 8);
return count;
}
int wordEigenvalues(DrawingContext *ctx, const char *word, double *out)
{
gsl_matrix *tmp = getTempMatrix(ctx->ws);
gsl_vector *ev = getTempVector(ctx->ws);
gsl_matrix **gen = getTempMatrices(ctx->ws, 6);
initializeTriangleGeneratorsCurrent(gen, ctx);
gsl_matrix_set_identity(tmp);
for(int i = 0; i < strlen(word); i++) {
if(word[i] == ' ')
continue;
multiply_right(tmp, gen[word[i]-'a'], ctx->ws);
}
int count = real_eigenvalues(tmp, ev, ctx->ws);
LOOP(i) out[i] = gsl_vector_get(ev, i);
releaseTempMatrices(ctx->ws, 7);
releaseTempVectors(ctx->ws, 1);
return count;
}
// level 1: the elementary drawing functions, drawPoint, drawSegment2d
void drawPoint(DrawingContext *ctx, point_t p)
{
cairo_t *C = ctx->cairo;
cairo_save(C);
cairo_arc(C, p.x, p.y, 3.0/ctx->dim->scalefactor, 0, 2*M_PI);
cairo_close_path(C);
cairo_fill(C);
cairo_restore(C);
/*
cairo_save(C);
cairo_move_to(C, p.x, p.y);
cairo_close_path(C);
cairo_set_line_cap(C, CAIRO_LINE_CAP_ROUND);
cairo_set_line_width(C, 10.0/ctx->dim->scalefactor);
cairo_stroke(C);
cairo_restore(C);
*/
}
void drawSegment2d(DrawingContext *ctx, point_t a, point_t b)
{
cairo_t *C = ctx->cairo;
cairo_move_to(C, a.x, a.y);
cairo_line_to(C, b.x, b.y);
cairo_stroke(C);
}
// level 2: drawVector, drawCovector, drawSegment
point_t vectorToPoint(DrawingContext *ctx, vector_t in)
{
double x[3];
point_t out;
LOOP(i) x[i] = 0.0;
LOOP(i) LOOP(j) x[i] += gsl_matrix_get(ctx->cob, i, j) * in.x[j];
out.x = x[0] / x[2];
out.y = x[1] / x[2];
return out;
}
void drawVector(DrawingContext *ctx, vector_t v)
{
drawPoint(ctx, vectorToPoint(ctx, v));
}
static void drawImplicitLine(DrawingContext *ctx, double a, double b, double c)
{
double norm, lambda;
point_t m, s, xminus, xplus;
m.x = ctx->dim->center_x;
m.y = ctx->dim->center_y;
lambda = (a*m.x + b*m.y + c)/(a*a + b*b);
s.x = m.x - lambda*a;
s.y = m.y - lambda*b;
norm = sqrt(a*a + b*b);
xminus.x = s.x - ctx->dim->radius * b / norm;
xminus.y = s.y + ctx->dim->radius * a / norm;
xplus.x = s.x + ctx->dim->radius * b / norm;
xplus.y = s.y - ctx->dim->radius * a / norm;
drawSegment2d(ctx, xminus, xplus);
}
void drawCovector(DrawingContext *ctx, vector_t v)
{
double x[3];
double cofactor;
LOOP(i) x[i] = 0.0;
LOOP(i) LOOP(j) {
cofactor = gsl_matrix_get(ctx->cob, (i+1)%3, (j+1)%3) * gsl_matrix_get(ctx->cob, (i+2)%3, (j+2)%3)
- gsl_matrix_get(ctx->cob, (i+1)%3, (j+2)%3) * gsl_matrix_get(ctx->cob, (i+2)%3, (j+1)%3);
x[i] += cofactor * v.x[j];
}
drawImplicitLine(ctx, x[0], x[1], x[2]);
}
void drawSegment(DrawingContext *ctx, vector_t a, vector_t b)
{
drawSegment2d(ctx, vectorToPoint(ctx, a), vectorToPoint(ctx, b));
}
void drawSegmentWith(DrawingContext *ctx, vector_t a, vector_t b, vector_t line, int contains)
{
point_t a_ = vectorToPoint(ctx,a);
point_t b_ = vectorToPoint(ctx,b);
double x[3];
double cofactor;
double r, tline, tminus, tplus;
double coeff0, coeff1, coeff2;
point_t m, xminus, xplus;
// multiply line with inverse of cob to get it as implicit line x in chart
LOOP(i) x[i] = 0.0;
LOOP(i) LOOP(j) {
cofactor = gsl_matrix_get(ctx->cob, (i+1)%3, (j+1)%3) * gsl_matrix_get(ctx->cob, (i+2)%3, (j+2)%3)
- gsl_matrix_get(ctx->cob, (i+1)%3, (j+2)%3) * gsl_matrix_get(ctx->cob, (i+2)%3, (j+1)%3);
x[i] += cofactor * line.x[j];
}
// t = parameter on segment of intersection with line, s(t) = a + (b-a)*t
tline = -(a_.x*x[0] + a_.y*x[1] + x[2])/((b_.x - a_.x)*x[0] + (b_.y - a_.y)*x[1]);
/*
printf("tline: %f\n", tline);
point_t s;
s.x = a_.x - (b_.x-a_.x)*tline;
s.y = a_.y - (b_.y-a_.y)*tline;
drawPoint(ctx, s);
*/
if((tline < 0 || tline > 1) != contains) {
drawSegment2d(ctx, a_, b_);
} else { // need to draw complementary segment
// find t so that s(t) is at radius r from center, |s(t)-m| = r
m.x = ctx->dim->center_x;
m.y = ctx->dim->center_y;
r = ctx->dim->radius;
// equation is coeff2 t^2 + 2 coeff1 t + coeff0 = 0
coeff0 = (a_.x - m.x)*(a_.x - m.x) + (a_.y - m.y)*(a_.y - m.y) - r*r;
coeff1 = (a_.x - m.x)*(b_.x - a_.x) + (a_.y - m.y)*(b_.y - a_.y);
coeff2 = (b_.x - a_.x)*(b_.x - a_.x) + (b_.y - a_.y)*(b_.y - a_.y);
if(coeff1*coeff1 - coeff0*coeff2 <= 0)
return;
tplus = (- coeff1 + sqrt(coeff1*coeff1 - coeff0*coeff2))/coeff2;
tminus = (- coeff1 - sqrt(coeff1*coeff1 - coeff0*coeff2))/coeff2;
xplus.x = a_.x + tplus * (b_.x - a_.x);
xplus.y = a_.y + tplus * (b_.y - a_.y);
xminus.x = a_.x + tminus * (b_.x - a_.x);
xminus.y = a_.y + tminus * (b_.y - a_.y);
if(tplus > 1)
drawSegment2d(ctx, b_, xplus);
if(tminus < 0)
drawSegment2d(ctx, a_, xminus);
}
}
// level 3: boxes and polygons
void drawPolygon(DrawingContext *ctx, int segments, int sides, ...)
{
va_list args;
vector_t first, prev, current;
va_start(args, sides);
first = va_arg(args, vector_t);
current = first;
for(int i = 0; i < sides-1; i++) {
prev = current;
current = va_arg(args, vector_t);
if(segments)
drawSegment(ctx, prev, current);
else
drawCovector(ctx, cross(prev, current));
}
if(segments)
drawSegment(ctx, current, first);
else
drawCovector(ctx, cross(current, first));
va_end(args);
}
void drawTriangle(DrawingContext *ctx, const char *word)
{
vector_t p[3];
fixedPoints(ctx, word, p);
drawPolygon(ctx, 1, 3, p[0], p[1], p[2]);
}
void drawBox(DrawingContext *ctx, const char *word1, const char *word2)
{
vector_t p[2][3],i[2];
fixedPoints(ctx, word1, p[0]);
fixedPoints(ctx, word2, p[1]);
// intersect attracting line with neutral line of the other element
for(int j = 0; j < 2; j++)
i[j] = cross(cross(p[j%2][0],p[j%2][1]),cross(p[(j+1)%2][0],p[(j+1)%2][2]));
drawPolygon(ctx, 1, 4, p[0][0], i[0], p[1][0], i[1]);
}
void drawBoxLines(DrawingContext *ctx, const char *word1, const char *word2)
{
vector_t p[2][3],i[2];
fixedPoints(ctx, word1, p[0]);
fixedPoints(ctx, word2, p[1]);
// intersect attracting line with neutral line of the other element
for(int j = 0; j < 2; j++)
i[j] = cross(cross(p[j%2][0],p[j%2][1]),cross(p[(j+1)%2][0],p[(j+1)%2][2]));
drawPolygon(ctx, 0, 4, p[0][0], i[0], p[1][0], i[1]);
}
void drawBoxStd(DrawingContext *ctx, const char *word, char base)
{
char word1[100];
char word2[100];
int len = strlen(word);
if(len*2 + 4 > 100)
return;
for(int i = 0; i < len; i++) {
word1[i] = word1[2*len+2-i] = word[i];
word2[i] = word2[2*len+2-i] = word[i];
}
word1[2*len+3] = 0;
word2[2*len+3] = 0;
LOOP(i) word1[len+i] = (base-'A'+6+i+1)%3+'a';
LOOP(i) word2[len+i] = (base-'A'+6-i-1)%3+'a';
// printf("Words: %s %s\n", word1, word2);
drawBox(ctx, word1, word2);
}
void drawRotationOrbitFrame(DrawingContext *ctx, gsl_matrix *frame, vector_t start)
{
vector_t v[3], w;
point_t p;
double parameter, startangle;
int iterations = 2000;
gsl_matrix *inverse = getTempMatrix(ctx->ws);
gsl_vector *start_v = getTempVector(ctx->ws);
gsl_vector *start_in_frame = getTempVector(ctx->ws);
cairo_t *C = ctx->cairo;
// computeRotationMatrixFrame(ctx, frame, word);
LOOP(i) LOOP(j) v[i].x[j] = gsl_matrix_get(frame, j, i);
LOOP(i) gsl_vector_set(start_v, i, start.x[i]);
solve(frame, start_v, start_in_frame, ctx->ws);
parameter = gsl_vector_get(start_in_frame, 2);
parameter /= sqrt(gsl_vector_get(start_in_frame, 0)*gsl_vector_get(start_in_frame, 0) +
gsl_vector_get(start_in_frame, 1)*gsl_vector_get(start_in_frame, 1));
startangle = atan2(gsl_vector_get(start_in_frame, 1), gsl_vector_get(start_in_frame, 0));
int previous_inside = 0;
for(int k = 0; k <= iterations; k++) {
LOOP(i) w.x[i] = parameter * v[2].x[i] + cos(2*k*M_PI/iterations) * v[0].x[i] + sin(2*k*M_PI/iterations) * v[1].x[i];
p = vectorToPoint(ctx, w);
if(isInsideBB(ctx, p)) {
if(!previous_inside)
cairo_move_to(C, p.x, p.y);
else
cairo_line_to(C, p.x, p.y);
previous_inside = 1;
} else {
previous_inside = 0;
}
}
cairo_stroke(C);
releaseTempMatrices(ctx->ws, 1);
releaseTempVectors(ctx->ws, 2);
}
void drawRotationOrbit(DrawingContext *ctx, const char *word, vector_t start)
{
gsl_matrix *frame = getTempMatrix(ctx->ws);
computeRotationMatrixFrame(ctx, frame, word);
drawRotationOrbitFrame(ctx, frame, start);
releaseTempMatrices(ctx->ws, 1);
}
void drawDualRotationOrbit(DrawingContext *ctx, const char *word, vector_t start)
{
vector_t v[3], w;
point_t p;
double parameter, startangle;
int iterations = 200;
gsl_matrix *frame = getTempMatrix(ctx->ws);
gsl_matrix *inverse = getTempMatrix(ctx->ws);
gsl_vector *start_v = getTempVector(ctx->ws);
gsl_vector *start_in_frame = getTempVector(ctx->ws);
cairo_t *C = ctx->cairo;
computeRotationMatrixFrame(ctx, frame, word);
LOOP(i) LOOP(j) v[i].x[j] = gsl_matrix_get(frame, j, i);
LOOP(i) gsl_vector_set(start_v, i, start.x[i]);
// solve(frame, start_v, start_in_frame, ctx->ws);
gsl_blas_dgemv(CblasTrans, 1.0, frame, start_v, 0.0, start_in_frame);
parameter = sqrt(gsl_vector_get(start_in_frame, 0)*gsl_vector_get(start_in_frame, 0) +
gsl_vector_get(start_in_frame, 1)*gsl_vector_get(start_in_frame, 1));
parameter /= gsl_vector_get(start_in_frame, 2);
startangle = atan2(gsl_vector_get(start_in_frame, 1), gsl_vector_get(start_in_frame, 0));
int previous_inside = 0;
for(int k = 0; k <= iterations; k++) {
LOOP(i) w.x[i] = parameter * v[2].x[i] + cos(2*k*M_PI/iterations) * v[0].x[i] + sin(2*k*M_PI/iterations) * v[1].x[i];
p = vectorToPoint(ctx, w);
if(isInsideBB(ctx, p)) {
if(!previous_inside)
cairo_move_to(C, p.x, p.y);
else
cairo_line_to(C, p.x, p.y);
previous_inside = 1;
} else {
previous_inside = 0;
}
}
cairo_stroke(C);
releaseTempMatrices(ctx->ws, 2);
releaseTempVectors(ctx->ws, 2);
}
void drawArcWithOutput(DrawingContext *ctx, const char *word, vector_t start, vector_type_t starttype, vector_t end, vector_type_t endtype, vector_t third, int contain, vector_t *start_vector_out, vector_t *end_vector_out, int dontdraw)
{
vector_t v[3], w, w_;
point_t p, p_;
double radius, angle_start, angle_end, angle_third, angle, angle_end_delta, sign, angle_start_final, angle_end_final, angle_end_other;
int iterations = 200;
gsl_matrix *frame = getTempMatrix(ctx->ws);
gsl_matrix *inverse = getTempMatrix(ctx->ws);
gsl_vector *vector = getTempVector(ctx->ws);
gsl_vector *vector_in_frame = getTempVector(ctx->ws);
cairo_t *C = ctx->cairo;
computeRotationMatrixFrame(ctx, frame, word);
LOOP(i) LOOP(j) v[i].x[j] = gsl_matrix_get(frame, j, i);
LOOP(i) gsl_vector_set(vector, i, start.x[i]);
if(starttype == VT_POINT) {
solve(frame, vector, vector_in_frame, ctx->ws);
radius = sqrt(gsl_vector_get(vector_in_frame, 0)*gsl_vector_get(vector_in_frame, 0) +
gsl_vector_get(vector_in_frame, 1)*gsl_vector_get(vector_in_frame, 1));
radius /= fabs(gsl_vector_get(vector_in_frame, 2));
angle_start = atan2(gsl_vector_get(vector_in_frame, 1)/gsl_vector_get(vector_in_frame, 2),
gsl_vector_get(vector_in_frame, 0)/gsl_vector_get(vector_in_frame, 2));
} else {
gsl_blas_dgemv(CblasTrans, 1.0, frame, vector, 0.0, vector_in_frame);
radius = fabs(gsl_vector_get(vector_in_frame, 2));
radius /= sqrt(gsl_vector_get(vector_in_frame, 0)*gsl_vector_get(vector_in_frame, 0) +
gsl_vector_get(vector_in_frame, 1)*gsl_vector_get(vector_in_frame, 1));
angle_start = atan2(gsl_vector_get(vector_in_frame, 1)/gsl_vector_get(vector_in_frame, 2),
gsl_vector_get(vector_in_frame, 0)/gsl_vector_get(vector_in_frame, 2));
}
LOOP(i) gsl_vector_set(vector, i, third.x[i]);
solve(frame, vector, vector_in_frame, ctx->ws);
angle_third = atan2(gsl_vector_get(vector_in_frame, 1)/gsl_vector_get(vector_in_frame, 2),
gsl_vector_get(vector_in_frame, 0)/gsl_vector_get(vector_in_frame, 2));
LOOP(i) gsl_vector_set(vector, i, end.x[i]);
if(endtype == VT_POINT) {
solve(frame, vector, vector_in_frame, ctx->ws);
angle_end = atan2(gsl_vector_get(vector_in_frame, 1)/gsl_vector_get(vector_in_frame, 2),
gsl_vector_get(vector_in_frame, 0)/gsl_vector_get(vector_in_frame, 2));
} else {
gsl_blas_dgemv(CblasTrans, 1.0, frame, vector, 0.0, vector_in_frame);
// this is only the average angle
angle_end = atan2(gsl_vector_get(vector_in_frame, 1)/gsl_vector_get(vector_in_frame, 2),
gsl_vector_get(vector_in_frame, 0)/gsl_vector_get(vector_in_frame, 2));
angle_end_delta = acos(-fabs(gsl_vector_get(vector_in_frame, 2))/radius/
sqrt(gsl_vector_get(vector_in_frame, 0)*gsl_vector_get(vector_in_frame, 0) +
gsl_vector_get(vector_in_frame, 1)*gsl_vector_get(vector_in_frame, 1)));
}
int previous_inside = 0;
for(int i = 0; i < 4; i++) {
angle_start_final = angle_start;
if(endtype == VT_POINT) {
angle_end_final = angle_end;
} else {
if(i >= 2) {
angle_end_final = angle_end - angle_end_delta;
angle_end_other = angle_end + angle_end_delta;
} else {
angle_end_final = angle_end + angle_end_delta;
angle_end_other = angle_end - angle_end_delta;
}
}
if(i%2)
FLIP(angle_start_final, angle_end_final);
if(endtype == VT_LINE && ANGLE_IN_INTERVAL(angle_start_final, angle_end_final, angle_end_other))
continue;
if(contain && !ANGLE_IN_INTERVAL(angle_start_final, angle_end_final, angle_third))
continue;
if(!contain && ANGLE_IN_INTERVAL(angle_start_final, angle_end_final, angle_third))
continue;
break;
}
// output the start end end point
/*
LOOP(i) w.x[i] = v[2].x[i] / radius + cos(angle_start_final) * v[0].x[i] + sin(angle_start_final) * v[1].x[i];
p = vectorToPoint(ctx, w);
LOOP(i) w.x[i] = v[2].x[i] / radius + cos(angle_end_final) * v[0].x[i] + sin(angle_end_final) * v[1].x[i];
p_ = vectorToPoint(ctx, w);
printf("\\draw (%f,%f) -- (%f,%f);\n", p.x, p.y, p_.x, p_.y);
*/
if(!dontdraw) {
for(int k = 0; k <= iterations; k++) {
angle = angle_start_final + (double)k/(double)iterations * ANGLE_DIFF(angle_end_final, angle_start_final);
LOOP(i) w.x[i] = v[2].x[i] / radius + cos(angle) * v[0].x[i] + sin(angle) * v[1].x[i];
p = vectorToPoint(ctx, w);
if(isInsideBB(ctx, p)) {
if(!previous_inside)
cairo_move_to(C, p.x, p.y);
else
cairo_line_to(C, p.x, p.y);
previous_inside = 1;
} else {
previous_inside = 0;
}
}
}
if(start_vector_out)
LOOP(i) start_vector_out->x[i] = v[2].x[i] / radius +
cos(angle_start_final) * v[0].x[i] +
sin(angle_start_final) * v[1].x[i];
if(end_vector_out)
LOOP(i) end_vector_out->x[i] = v[2].x[i] / radius +
cos(angle_end_final) * v[0].x[i] +
sin(angle_end_final) * v[1].x[i];
cairo_stroke(C);
releaseTempMatrices(ctx->ws, 2);
releaseTempVectors(ctx->ws, 2);
}
void drawArc(DrawingContext *ctx, const char *word, vector_t start, vector_type_t starttype, vector_t end, vector_type_t endtype, vector_t third, int contain)
{
drawArcWithOutput(ctx, word, start, starttype, end, endtype, third, contain, 0, 0, 0);
}
// level 4: draw the actual image components
void drawReflectors(DrawingContext *ctx)
{
vector_t v[3];
cairo_set_source_rgb(ctx->cairo, 0, 0, 0);
LOOP(i) LOOP(j) { v[i].x[j] = (i==j) ? 1.0 : 0.0; }
LOOP(i) drawVector(ctx, v[i]);
LOOP(i) LOOP(j) v[i].x[j] = gsl_matrix_get(ctx->cartan, i, j);
LOOP(i) drawCovector(ctx, v[i]);
}
void drawAttractors(DrawingContext *ctx)
{
int n = 3;
vector_t p[6][3];
vector_t l[6][3];
fixedPoints(ctx, "cba", p[0]);
fixedPoints(ctx, "bac", p[1]);
fixedPoints(ctx, "acb", p[2]);
fixedPoints(ctx, "a cab a", p[3]);
fixedPoints(ctx, "b abc b", p[4]);
fixedPoints(ctx, "c bca c", p[5]);
double color[6][3] = {{1,0,0},{0,0.7,0},{0,0,1},{0,1,1},{0,1,1},{0,1,1}};
for(int i = 0; i < n; i++) LOOP(j) l[i][j] = cross(p[i][(3-j)%3], p[i][(4-j)%3]);
for(int i = 0; i < n; i++) LOOP(j) {
cairo_set_source_rgb(ctx->cairo, color[i][0], color[i][1], color[i][2]);
drawVector(ctx, p[i][j]);
}
for(int i = 0; i < n; i++) LOOP(j) {
cairo_set_source_rgb(ctx->cairo, color[i][0], color[i][1], color[i][2]);
drawCovector(ctx, l[i][j]);
}
}
char *conjugate_word(const char *word, int modifier, const char *conj, char *buffer)
{
int wordlen = strlen(word);
int conjlen = strlen(conj);
for(int i = 0; i < conjlen; i++) {
buffer[i] = conj[i];
buffer[2*conjlen+wordlen-1-i] = conj[i];
}
for(int i = 0; i < wordlen; i++) {
if(word[i] == ' ')
buffer[conjlen+i] = word[i];
else
buffer[conjlen+i] = (word[i]+modifier-'a')%3 + 'a';
}
buffer[2*conjlen + wordlen] = 0;
return buffer;
}
void drawCurvedBox(DrawingContext *ctx, int base, const char *conj, int style)
{
vector_t p[11][3];
vector_t l[2][3];
vector_t corner1, corner2;
vector_t tmp1, tmp2;
char word[100];
int modifier = base - 'A';
conjugate_word("abc", modifier, conj, word);
fixedPoints(ctx, word, p[0]);
conjugate_word("bca", modifier, conj, word);
fixedPoints(ctx, word, p[1]);
conjugate_word("b abc b", modifier, conj, word);
fixedPoints(ctx, word, p[2]);
conjugate_word("ab abc ba", modifier, conj, word);
fixedPoints(ctx, word, p[3]);
conjugate_word("baca cab acab", modifier, conj, word);
fixedPoints(ctx, word, p[4]);
conjugate_word("abaca cab acaba", modifier, conj, word);
fixedPoints(ctx, word, p[5]);
conjugate_word("bca cab acb", modifier, conj, word);
fixedPoints(ctx, word, p[6]);
conjugate_word("abca cab acba", modifier, conj, word);
fixedPoints(ctx, word, p[7]);
conjugate_word("abacababa", modifier, conj, word);
fixedPoints(ctx, word, p[8]);
conjugate_word("bacabab", modifier, conj, word);
fixedPoints(ctx, word, p[9]);
conjugate_word("cab", modifier, conj, word);
fixedPoints(ctx, word, p[10]);
// conjugate_word("bca b abc b acb", modifier, conj, word);
// fixedPoints(ctx, word, p[6]);
// conjugate_word("bca baca cab acab acb", modifier, conj, word);
// fixedPoints(ctx, word, p[7]);
LOOP(j) l[0][j] = cross(p[0][(3-j)%3], p[0][(4-j)%3]);
LOOP(j) l[1][j] = cross(p[1][(3-j)%3], p[1][(4-j)%3]);
LOOP(j) l[10][j] = cross(p[10][(3-j)%3], p[10][(4-j)%3]);
// main conic
conjugate_word("ab", modifier, conj, word);
drawArcWithOutput(ctx, word, p[0][0], VT_POINT, p[2][0], VT_POINT, p[1][0], 0, &tmp1, &tmp2, style != 1);
if(style == 2)
drawSegmentWith(ctx, tmp1, tmp2, l[10][0], 0);
if(style == 3) {
drawVector(ctx, tmp1);
drawVector(ctx, tmp2);
}
conjugate_word("ab", modifier, conj, word);
drawArcWithOutput(ctx, word, p[1][0], VT_POINT, p[3][0], VT_POINT, p[0][0], 0, &tmp1, &tmp2, style != 1);
if(style == 2)
drawSegmentWith(ctx, tmp1, tmp2, l[10][0], 0);
if(style == 3) {
drawVector(ctx, tmp1);
drawVector(ctx, tmp2);
}
conjugate_word("bcabcb", modifier, conj, word);
drawArcWithOutput(ctx, word, p[2][0], VT_POINT, l[1][0], VT_LINE, p[6][0], 1, &tmp1, &tmp2, style != 1); // only 1st cutoff
if(style == 2)
drawSegmentWith(ctx, tmp1, tmp2, l[10][0], 0);
corner1 = tmp2;
conjugate_word("abcabcba", modifier, conj, word);
drawArcWithOutput(ctx, word, p[3][0], VT_POINT, l[0][0], VT_LINE, p[7][0], 1, &tmp1, &tmp2, style != 1); // only 1st cutoff
if(style == 2)
drawSegmentWith(ctx, tmp1, tmp2, l[10][0], 0);
corner2 = tmp2;
if(style == 1 || style == 2) {
drawSegmentWith(ctx, p[0][0], corner2, l[1][1], 0);
drawSegmentWith(ctx, p[1][0], corner1, l[0][1], 0);
} else if(style == 3)
{
drawVector(ctx, corner1);
drawVector(ctx, corner2);
}
}
groupelement_t *left(const char *word, groupelement_t *g)
{
int n = strlen(word);
for(int i = n-1; i >= 0; i--) {
if(word[i] == ' ')
continue;
g = g->adj[word[i]-'a'];
if(!g)
break;
}
return g;
}
void drawBoxes(DrawingContext *ctx)
{
gsl_matrix *tmp = getTempMatrix(ctx->ws);
gsl_matrix **gen = getTempMatrices(ctx->ws, 6);
gsl_matrix *frame = getTempMatrix(ctx->ws);
gsl_matrix *frame2 = getTempMatrix(ctx->ws);
gsl_vector *startpoint_drawbasis = getTempVector(ctx->ws);
gsl_vector *startpoint_globalbasis = getTempVector(ctx->ws);
gsl_matrix **elements = getTempMatrices(ctx->ws, ctx->n_group_elements);
cairo_t *C = ctx->cairo;
cairo_save(C);
cairo_set_line_width(C, 2.0/ctx->dim->scalefactor);
cairo_set_source_rgb(C, 0.6, 0.6, 0.6);
vector_t p[22][3];
vector_t fp[3];
vector_t l[22][3];
vector_t alpha[6];
vector_t ptmp[3];
vector_t start;
vector_t start2;
char word[100], word2[100];
fixedPoints(ctx, "cba", p[0]);
fixedPoints(ctx, "acb", p[1]);
fixedPoints(ctx, "bac", p[2]);
initializeTriangleGeneratorsCurrent(gen, ctx);
gsl_matrix_set_identity(elements[0]);
for(int i = 1; i < ctx->n_group_elements; i++) {
if(ctx->group[i].length % 2)
continue;
int letter = ROTATION_LETTER(ctx->group[i].letter, ctx->group[i].parent->letter);
multiply(gen[letter], elements[ctx->group[i].parent->parent->id], elements[i]);
}
gsl_vector_set(startpoint_drawbasis, 0, ctx->marking.x);
gsl_vector_set(startpoint_drawbasis, 1, ctx->marking.y);
gsl_vector_set(startpoint_drawbasis, 2, 1);
solve(ctx->cob, startpoint_drawbasis, startpoint_globalbasis, ctx->ws);
LOOP(i) start.x[i] = gsl_vector_get(startpoint_globalbasis, i);
/*
fixedPoints(ctx, "ABAB", fp);
drawCovector(ctx, cross(fp[0], fp[2]));
drawVector(ctx, fp[1]);
computeRotationMatrixFrame(ctx, frame, "c");
drawRotationOrbit(ctx, "c", start); // ba cb ab ac = C A c B
drawVector(ctx, start);
computeMatrix(ctx, tmp, "ABABc");
start2 = apply(tmp, start);
computeRotationMatrixFrame(ctx, frame, "c");
multiply(tmp, frame, frame2);
drawVector(ctx, start2);
drawRotationOrbitFrame(ctx, frame2, apply(tmp, start));
computeMatrix(ctx, tmp, "ABABcABABC");
start2 = apply(tmp, start);
computeRotationMatrixFrame(ctx, frame, "c");
multiply(tmp, frame, frame2);
drawVector(ctx, start2);
drawRotationOrbitFrame(ctx, frame2, apply(tmp, start));
*/
for(groupelement_t *cur = &ctx->group[180]; cur->parent; cur = cur->parent)
fputc('a'+cur->letter, stdout);
fputc('\n', stdout);
queue_t queue;
queue_init(&queue);
queue_put(&queue, 0);
int current;
groupelement_t *cur, *next;
for(int i = 0; i < ctx->n_group_elements_combinatorial; i++)
ctx->group[i].visited = 0;
cur = &ctx->group[0];
computeRotationMatrixFrame(ctx, frame, "c");
for(int i = 0; i < 1000; i++) {
if(ctx->group[i].length % 2 != 0)
continue;
multiply(elements[i], frame, frame2);
drawRotationOrbitFrame(ctx, frame2, apply(elements[i], start));
}
/*
while((current = queue_get(&queue)) != -1) {
cur = &ctx->group[current];
if(cur->visited > 4)
continue;
if(cur->id < ctx->n_group_elements) {
multiply(elements[cur->id], frame, frame2);
drawRotationOrbitFrame(ctx, frame2, apply(elements[cur->id], start));
}
next = left("ab ab", cur);
if(next && next->visited == 0) {
queue_put(&queue, next->id);
next->visited = cur->visited+1;
}
next = left("cbac cbac", cur);
if(next && next->visited == 0) {
queue_put(&queue, next->id);
next->visited = cur->visited+1;
}
next = left("cacbca cacbca", cur);
if(next && next->visited == 0) {
queue_put(&queue, next->id);
next->visited = cur->visited+1;
}
next = left("cabcacbcac cabcacbcac", cur);
if(next && next->visited == 0) {
queue_put(&queue, next->id);
next->visited = cur->visited+1;
}
next = left("acbcacba acbcacba", cur);
if(next && next->visited == 0) {
queue_put(&queue, next->id);
next->visited = cur->visited+1;
}
next = left("bcacbc bcacbc", cur);
if(next && next->visited == 0) {
queue_put(&queue, next->id);
next->visited = cur->visited+1;
}
}
*/
// drawRotationOrbit(ctx, "a", p[1][0]);
// drawRotationOrbit(ctx, "b", p[2][0]);
cairo_restore(C);
releaseTempMatrices(ctx->ws, 9 + ctx->n_group_elements);
releaseTempVectors(ctx->ws, 2);
}
void drawBoxes2(DrawingContext *ctx)
{
gsl_matrix *rot = getTempMatrix(ctx->ws);
gsl_matrix **gen = getTempMatrices(ctx->ws, 6);
cairo_t *C = ctx->cairo;
cairo_save(C);
initializeTriangleGeneratorsCurrent(gen, ctx);
vector_t p[4][3];
fixedPoints(ctx, "cba", p[0]);
fixedPoints(ctx, "acb", p[1]);
fixedPoints(ctx, "bac", p[2]);
cairo_set_line_width(C, 2.5/ctx->dim->scalefactor);
/*
cairo_set_source_rgb(C, 0, 1, 1);
fixedPoints(ctx, "abc", p[0]);
drawRotationOrbit(ctx, "bc", p[0][0]); // (r^{-1}a)^2 C / cbr
fixedPoints(ctx, "bca", p[0]);
drawRotationOrbit(ctx, "abca", p[0][0]); // (r^{-1}a)^2 C / cbr
cairo_set_source_rgb(C, 1, 1, 0);
fixedPoints(ctx, "abc", p[0]);
drawRotationOrbit(ctx, "cbcabc", p[0][0]); // (r^{-1}a)^4 C / cbac r-
fixedPoints(ctx, "bca", p[0]);
drawRotationOrbit(ctx, "acbcabca", p[0][0]); // (r^{-1}a)^4 C / cbac r-
cairo_set_source_rgb(C, 1, 0, 1);
fixedPoints(ctx, "abc", p[0]);
drawRotationOrbit(ctx, "cbacabcabc", p[0][0]); // (r^{-1}a)^6 C / cbacba
fixedPoints(ctx, "bca", p[0]);
drawRotationOrbit(ctx, "acbacabcabca", p[0][0]); // (r^{-1}a)^6 C / cbacba
*/
/*
cairo_set_source_rgb(C, 1, 0, 0);
fixedPoints(ctx, "bca", p[3]);
drawRotationOrbit(ctx, "bcabcb", p[3][0]); // cb c a bc
fixedPoints(ctx, "abc", p[3]);
drawRotationOrbit(ctx, "abcabcba", p[3][0]); // cb c a bc
cairo_set_source_rgb(C, 0, 1, 0);
fixedPoints(ctx, "ababcba", p[3]);
drawRotationOrbit(ctx, "ababcababa", p[3][0]); // cb c a bc
fixedPoints(ctx, "abc", p[3]);
drawRotationOrbit(ctx, "abcabcabacba", p[3][0]); // cb c a bc
*/
// drawRotationOrbit(ctx, "bc", p[3][0]);
cairo_set_source_rgb(C, 0, 0, 0);
drawCurvedBox(ctx, 'A', "", 2);
cairo_set_source_rgb(C, 0, 0, 1);
// drawCurvedBox(ctx, 'C', "ab", 2);
// drawCurvedBox(ctx, 'C', "ab", 2);
// drawCurvedBox(ctx, 'B', "bca", 2);
// drawCurvedBox(ctx, 'B', "ba", 2);
// drawCurvedBox(ctx, 'B', "bca", 2);
// drawCurvedBox(ctx, 'B', "abca", 2);
// drawCurvedBox(ctx, 'C', "abab", 2);
// drawCurvedBox(ctx, 'C', "ababab", 2);
if(ctx->mode >= 3 && ctx->mode != 4) {
cairo_set_source_rgb(C, 0, 0, 0);
drawCurvedBox(ctx, 'B', "", 2);
drawCurvedBox(ctx, 'C', "", 2);
}
if(ctx->mode == 4 || ctx->mode == 5) {
cairo_set_source_rgb(C, 0, 0.8, 0.5);
drawCurvedBox(ctx, 'A', "", 2);
drawCurvedBox(ctx, 'A', "bc", 2);
drawCurvedBox(ctx, 'A', "bcbc", 2);
drawCurvedBox(ctx, 'A', "bcbcbc", 2);
drawCurvedBox(ctx, 'A', "bcbcbcbc", 2);
drawCurvedBox(ctx, 'A', "bcbcbcbcbc", 2);
drawCurvedBox(ctx, 'A', "bcbcbcbcbcbc", 2);
// drawCurvedBox(ctx, 'A', "b", 2);
// drawCurvedBox(ctx, 'A', "bcb", 2);
// drawCurvedBox(ctx, 'A', "bcbcb", 2);
// drawCurvedBox(ctx, 'A', "bcbcbcb", 2);
// drawCurvedBox(ctx, 'A', "bcbcbcbcb", 2);
}
if(ctx->mode == 6) {
cairo_set_source_rgb(C, 0, 0.8, 0.5);
drawCurvedBox(ctx, 'C', "abababab", 2);
drawCurvedBox(ctx, 'C', "ababababab", 2);
}
if(ctx->mode >= 6) {
cairo_set_source_rgb(C, 0, 0.8, 0.5);
drawCurvedBox(ctx, 'C', "ab", 2);
drawCurvedBox(ctx, 'C', "abab", 2);
drawCurvedBox(ctx, 'C', "ababab", 2);
}
if(ctx->mode >= 8) {
cairo_set_source_rgb(C, 0, 0.8, 0.5);
drawCurvedBox(ctx, 'C', "b", 2);
drawCurvedBox(ctx, 'C', "bab", 2);
drawCurvedBox(ctx, 'C', "babab", 2);
}
if(ctx->mode >= 9) {
cairo_set_source_rgb(C, 0.8, 0, 0.5);
drawCurvedBox(ctx, 'B', "abca", 2);
drawCurvedBox(ctx, 'B', "abcaca", 2);
drawCurvedBox(ctx, 'B', "abcacaca", 2);
drawCurvedBox(ctx, 'B', "aba", 2);
drawCurvedBox(ctx, 'B', "abaca", 2);
drawCurvedBox(ctx, 'B', "abacaca", 2);
drawCurvedBox(ctx, 'B', "ababca", 2);
drawCurvedBox(ctx, 'B', "ababcaca", 2);
drawCurvedBox(ctx, 'B', "ababcacaca", 2);
drawCurvedBox(ctx, 'B', "ababa", 2);
drawCurvedBox(ctx, 'B', "ababaca", 2);
drawCurvedBox(ctx, 'B', "ababacaca", 2);
drawCurvedBox(ctx, 'B', "abababca", 2);
drawCurvedBox(ctx, 'B', "abababcaca", 2);
drawCurvedBox(ctx, 'B', "abababcacaca", 2);
drawCurvedBox(ctx, 'B', "abababa", 2);
drawCurvedBox(ctx, 'B', "abababaca", 2);
drawCurvedBox(ctx, 'B', "abababacaca", 2);
drawCurvedBox(ctx, 'B', "bca", 2);
drawCurvedBox(ctx, 'B', "bcaca", 2);
drawCurvedBox(ctx, 'B', "bcacaca", 2);
drawCurvedBox(ctx, 'B', "ba", 2);
drawCurvedBox(ctx, 'B', "baca", 2);
drawCurvedBox(ctx, 'B', "bacaca", 2);
drawCurvedBox(ctx, 'B', "babca", 2);
drawCurvedBox(ctx, 'B', "babcaca", 2);
drawCurvedBox(ctx, 'B', "babcacaca", 2);
drawCurvedBox(ctx, 'B', "baba", 2);
drawCurvedBox(ctx, 'B', "babaca", 2);
drawCurvedBox(ctx, 'B', "babacaca", 2);
drawCurvedBox(ctx, 'B', "bababca", 2);
drawCurvedBox(ctx, 'B', "bababcaca", 2);
drawCurvedBox(ctx, 'B', "bababcacaca", 2);
drawCurvedBox(ctx, 'B', "bababa", 2);
drawCurvedBox(ctx, 'B', "bababaca", 2);
drawCurvedBox(ctx, 'B', "bababacaca", 2);
}
// drawCurvedBox(ctx, 'B', "abab", 2);
// drawCurvedBox(ctx, 'B', "ababab", 2);
// drawCurvedBox(ctx, 'C', "b", 1);
// drawCurvedBox(ctx, 'C', "ababab");
// drawCurvedBox(ctx, 'C', "abababab");
// drawCurvedBox(ctx, 'C', "b");
// drawCurvedBox(ctx, 'C', "bab");
// drawCurvedBox(ctx, 'C', "babab");
// drawCurvedBox(ctx, 'C', "abababab");
// drawCurvedBox(ctx, 'C', "ababababab");
// drawCurvedBox(ctx, 'C', "abababababab");
cairo_set_source_rgb(C, 0, 0, 1);
// drawCurvedBox(ctx, 'B', "abca", 1);
// drawCurvedBox(ctx, 'B', "aba", 1);
// drawCurvedBox(ctx, 'B', "bca", 1);
// drawCurvedBox(ctx, 'B', "ba", 1);
// drawCurvedBox(ctx, 'B', "abaca");
// drawCurvedBox(ctx, 'B', "abcaca");
// drawCurvedBox(ctx, 'B', "abaca");
// drawCurvedBox(ctx, 'B', "aba");
// drawCurvedBox(ctx, 'B', "ababca");
// drawCurvedBox(ctx, 'B', "ababa");
// drawCurvedBox(ctx, 'B', "abababca");
// drawCurvedBox(ctx, 'B', "abababa");
// drawCurvedBox(ctx, 'B', "bca");
// drawCurvedBox(ctx, 'B', "ba");
// drawCurvedBox(ctx, 'B', "babca");
// drawCurvedBox(ctx, 'B', "baba");
// drawCurvedBox(ctx, 'B', "bababca");
// drawCurvedBox(ctx, 'B', "bababa");
cairo_set_source_rgb(C, 1, 0, 1);
// drawCurvedBox(ctx, 'A', "abcac");
// drawCurvedBox(ctx, 'A', "bcac");
// drawCurvedBox(ctx, 'A', "abcacbc");
// drawCurvedBox(ctx, 'A', "abcac");
// drawCurvedBox(ctx, 'A', "abcacbc");
// drawCurvedBox(ctx, 'A', "ababc");
// drawCurvedBox(ctx, 'A', "ababcbc");
/* drawCurvedBox(ctx, 'A', "abcac");
drawCurvedBox(ctx, 'A', "ababc");
drawCurvedBox(ctx, 'A', "abac");
drawCurvedBox(ctx, 'A', "ababcabc");
drawCurvedBox(ctx, 'A', "ababcac");
drawCurvedBox(ctx, 'A', "abababc");
drawCurvedBox(ctx, 'A', "ababac");
drawCurvedBox(ctx, 'A', "abababcabc");
drawCurvedBox(ctx, 'A', "abababcac");
drawCurvedBox(ctx, 'A', "ababababc");
drawCurvedBox(ctx, 'A', "abababac");
drawCurvedBox(ctx, 'A', "bcabc");
drawCurvedBox(ctx, 'A', "bcac");
drawCurvedBox(ctx, 'A', "babc");
drawCurvedBox(ctx, 'A', "bac");
drawCurvedBox(ctx, 'A', "babcabc");
drawCurvedBox(ctx, 'A', "babcac");
drawCurvedBox(ctx, 'A', "bababc");
drawCurvedBox(ctx, 'A', "babac");
drawCurvedBox(ctx, 'A', "bababcabc");
drawCurvedBox(ctx, 'A', "bababcac");
drawCurvedBox(ctx, 'A', "babababc");
drawCurvedBox(ctx, 'A', "bababac");
*/
// drawCurvedBox(ctx, 'A', "ababc");
cairo_set_source_rgb(C, 1, 0.5, 0);
// drawCurvedBox(ctx, 'C', "abcabcab");
// drawCurvedBox(ctx, 'B', "bca");
// drawCurvedBox(ctx, 'B', "abcacaca");
// drawCurvedBox(ctx, 'B', "abaca");
// drawCurvedBox(ctx, 'B', "aba");
// drawCurvedBox(ctx, 'B', "aba");
// drawCurvedBox(ctx, 'B', "abcaca");
// drawCurvedBox(ctx, 'B', "aba");
// drawCurvedBox(ctx, 'B', "abaca");
cairo_set_source_rgb(C, 0, 0, 1);
// drawCurvedBox(ctx, 'A', "abcabc");
// drawCurvedBox(ctx, 'A', "abcabcbc");
// drawCurvedBox(ctx, 'A', "abcacbc");
// drawCurvedBox(ctx, 'A', "abcac");
// drawCurvedBox(ctx, 'A', "bcabc");
// drawCurvedBox(ctx, 'A', "bcabcbc");
// drawCurvedBox(ctx, 'A', "bcacbc");
// drawCurvedBox(ctx, 'A', "bcac");
// drawCurvedBox(ctx, 'A', "abcac");
// drawCurvedBox(ctx, 'A', "abcac");
// drawCurvedBox(ctx, 'A', "abcacbc");
// drawCurvedBox(ctx, 'B', "abacaca");
// drawCurvedBox(ctx, 'B', "abacacaca");
cairo_set_source_rgb(C, 0, 1, 1);
fixedPoints(ctx, "ababcba", p[0]);
// drawRotationOrbit(ctx, "abcaba", p[0][0]);
/*
int p = 9;
vector_t fp[3][3],neutral_line[3],reflection_line[p],star[2*p],outer[2];
vector_t rotation_line = {0,0,1};
cairo_set_line_width(C, 1.5/ctx->dim->scalefactor);
cairo_set_source_rgb(C, 1, 0, 0);
multiply(gen[0], gen[1], rot);
fixedPoints(ctx, "abc", fp[0]);
fixedPoints(ctx, "bca", fp[1]);
fixedPoints(ctx, "cab", fp[2]);
LOOP(i) neutral_line[i] = cross(fp[i][0], fp[i][2]);
LOOP(j) reflection_line[0].x[j] = gsl_matrix_get(ctx->cartan, 0, j);
star[0] = cross(neutral_line[0],reflection_line[0]);
star[p] = cross(neutral_line[2],reflection_line[0]);
for(int j = 1; j < p; j++) {
reflection_line[j] = apply_transpose(rot, reflection_line[j-1]);
star[j] = apply(rot, star[j-1]);
star[j+p] = apply(rot, star[j+p-1]);
}
outer[0] = cross(neutral_line[0],reflection_line[5]);
outer[1] = cross(neutral_line[2],reflection_line[5]);
for(int j = 0; j < 8; j++) {
drawVector(ctx, star[j]);
drawVector(ctx, star[p+j]);
}
for(int j = 0; j < 7; j++) {
// if(j == 3)
// continue;
drawSegment(ctx, star[j%p], star[(j+1)%p+p]);
drawSegment(ctx, star[(j+1)%p+p], star[(j+1)%p]);
drawSegment(ctx, star[(j+1)%p], star[j%p+p]);
drawSegment(ctx, star[j%p+p], star[j%p]);
}
cairo_set_source_rgb(C, 0, 0, 1);
drawVector(ctx,outer[0]);
drawVector(ctx,outer[1]);
// drawCovector(ctx, neutral_line[0]);
// drawCovector(ctx, neutral_line[2]);
drawSegment(ctx, star[0], star[p]);
drawSegment(ctx, star[p], outer[1]);
drawSegment(ctx, outer[1], outer[0]);
drawSegment(ctx, outer[0], star[0]);
cairo_set_source_rgb(C, 0, 0, 0);
drawCovector(ctx, rotation_line);
*/
// for(int j = 0; j < 5; j++)
// drawCovector(ctx, reflection_line[j]);
// intersect attracting line with neutral line of the other element
// for(int j = 0; j < 2; j++)
// i[j] = cross(cross(p[j%2][0],p[j%2][1]),cross(p[(j+1)%2][0],p[(j+1)%2][2]));
// drawPolygon(ctx, 1, 4, p[0][0], i[0], p[1][0], i[1]);
// abc, ababcba, abababcbaba, ..., cab
// bca, acaba, abacababa, ..., babcb
/*
vector_t v[4][3];
vector_t i[4];
fixedPoints(ctx, "abc", v[0]);
fixedPoints(ctx, "bca", v[1]);
fixedPoints(ctx, "cab", v[2]);
fixedPoints(ctx, "acaba", v[3]);
i[0] = cross(cross(v[0][1],v[0][2]),cross(v[2][0],v[2][2]));
i[1] = cross(cross(v[1][1],v[1][2]),cross(v[3][0],v[3][2]));
i[2] = cross(cross(v[2][0],v[2][2]),cross(v[3][0],v[3][2]));
i[3] = cross(cross(v[0][0],v[0][2]),cross(v[1][0],v[1][2]));
// cairo_set_source_rgb(C, 0, 0, 1);
// drawPolygon(ctx, 1, 6, v[2][2], v[1][1], v[0][1], v[3][2], v[3][1], v[2][1]);
// drawPolygon(ctx, 1, 6, v[1][2], i[1], i[2], i[0], v[0][2], i[3]);
cairo_set_line_width(C, 1.5/ctx->dim->scalefactor);
cairo_set_source_rgb(C, 1, 0, 0);
// drawBox(ctx, "bca", "abc");
drawCurvedBox(ctx, 'A', "");
*/
/*
cairo_set_source_rgb(C, 1, 0.5, 0);
drawBox(ctx, "ab abc ba", "ab cab ba");
drawBox(ctx, "abab abc baba", "abab cab baba");
drawBox(ctx, "ababab abc bababa", "ababab cab bababa");
drawBox(ctx, "abababab abc babababa", "abababab cab babababa");
// drawBox(ctx, "ababababab abc bababababa", "ababababab cab bababababa");
drawBox(ctx, "b abc b", "b cab b");
drawBox(ctx, "bab abc bab", "bab cab bab");
drawBox(ctx, "babab abc babab", "babab cab babab");
drawBox(ctx, "bababab abc bababab", "bababab cab bababab");
// drawBox(ctx, "babababab abc babababab", "babababab cab babababab");
*/
cairo_restore(C);
releaseTempMatrices(ctx->ws, 7);
}
void drawRotatedReflectors(DrawingContext *ctx)
{
gsl_matrix *rot = getTempMatrix(ctx->ws);
gsl_matrix **gen = getTempMatrices(ctx->ws, 6);
cairo_t *C = ctx->cairo;
vector_t fp[3], fp2[3];
vector_t w;
vector_t v[3];
cairo_set_source_rgb(C, 0.7, 0.7, 0.7);
initializeTriangleGeneratorsCurrent(gen, ctx);
LOOP(i) LOOP(j) v[i].x[j] = gsl_matrix_get(ctx->cartan, i, j);
multiply(gen[0], gen[1], rot);
for(int j = 0; j < ctx->p[2]; j++) {
drawCovector(ctx, v[0]);
v[0] = apply_transpose(rot, v[0]);
}
LOOP(i) LOOP(j) { v[i].x[j] = (i==j) ? 1.0 : 0.0; }
for(int j = 0; j < ctx->p[2]; j++) {
drawVector(ctx, v[0]);
v[0] = apply(rot, v[0]);
}
fixedPoints(ctx, "cab", fp);
fixedPoints(ctx, "cacabac", fp2);
drawRotationOrbit(ctx, "ac", fp[0]);
releaseTempMatrices(ctx->ws, 7);
}
void drawDualLimitCurve(DrawingContext *ctx)
{
cairo_t *C = ctx->cairo;
cairo_save(C);
cairo_set_source_rgb(C, 0.5, 0.5, 1);
gsl_matrix **gen = getTempMatrices(ctx->ws, 6);
gsl_matrix **elements = getTempMatrices(ctx->ws, ctx->n_group_elements);
// wordEigenvalues(ctx, "abc", ev);
// LOOP(i) LOOP(j) gsl_matrix_set(coxeter_fixedpoints, j, i, cox[0][i].x[j]);
initializeTriangleGeneratorsCurrent(gen, ctx);
gsl_matrix_set_identity(elements[0]);
for(int i = 1; i < ctx->n_group_elements; i++)
multiply(gen[ctx->group[i].letter], elements[ctx->group[i].parent->id], elements[i]);
vector_t p[3], l[3], v;
fixedPoints(ctx, "cba", p);
drawVector(ctx, p[0]);
drawVector(ctx, p[1]);
drawVector(ctx, p[2]);
LOOP(i) l[i] = cross(p[(i+1)%3], p[(i+2)%3]);
for(int i = 0; i < ctx->n_group_elements; i++) {
v = apply_transpose(elements[i], l[0]);
drawCovector(ctx, v);
}
releaseTempMatrices(ctx->ws, 3 + ctx->n_group_elements);
// releaseTempVectors(ctx->ws, 4);
/*
fixedPoints(ctx, "ab abc ba", p[1]);
fixedPoints(ctx, "abab abc baba", p[2]);
fixedPoints(ctx, "ababab abc bababa", p[3]);
fixedPoints(ctx, "abababab abc babababa", p[4]);
fixedPoints(ctx, "babababa abc abababab", p[5]);
fixedPoints(ctx, "bababa abc ababab", p[6]);
fixedPoints(ctx, "baba abc abab", p[7]);
fixedPoints(ctx, "ba abc ab", p[8]);
fixedPoints(ctx, "bca", p[9]);
fixedPoints(ctx, "ab bca ba", p[10]);
fixedPoints(ctx, "abab bca baba", p[11]);
fixedPoints(ctx, "ababab bca bababa", p[12]);
fixedPoints(ctx, "abababab bca babababa", p[13]);
fixedPoints(ctx, "babababa bca abababab", p[14]);
fixedPoints(ctx, "bababa bca ababab", p[15]);
fixedPoints(ctx, "baba bca abab", p[16]);
fixedPoints(ctx, "ba bca ab", p[17]);
*/
/*
fixedPoints(ctx, "abc", p[0]);
fixedPoints(ctx, "ac abc ca", p[1]);
fixedPoints(ctx, "acac abc caca", p[2]);
fixedPoints(ctx, "acacac abc cacaca", p[3]);
fixedPoints(ctx, "acacacac abc cacacaca", p[4]);
fixedPoints(ctx, "cacacaca abc acacacac", p[5]);
fixedPoints(ctx, "cacaca abc acacac", p[6]);
fixedPoints(ctx, "caca abc acac", p[7]);
fixedPoints(ctx, "ca abc ac", p[8]);
fixedPoints(ctx, "bca", p[9]);
fixedPoints(ctx, "ac bca ca", p[10]);
fixedPoints(ctx, "acac bca caca", p[11]);
fixedPoints(ctx, "acacac bca cacaca", p[12]);
fixedPoints(ctx, "acacacac bca cacacaca", p[13]);
fixedPoints(ctx, "cacacaca bca acacacac", p[14]);
fixedPoints(ctx, "cacaca bca acacac", p[15]);
fixedPoints(ctx, "caca bca acac", p[16]);
fixedPoints(ctx, "ca bca ac", p[17]);
*/
/*
fixedPoints(ctx, "cab", p[2]);
fixedPoints(ctx, "b abc b", p[3]);
fixedPoints(ctx, "c bca c", p[4]);
fixedPoints(ctx, "a cab a", p[5]);
*/
/*
for(int i = 0; i < n; i++) {
LOOP(j) l[i][j] = cross(p[i][(3-j)%3], p[i][(4-j)%3]);
drawCovector(ctx, l[i][0]);
// drawCovector(ctx, l[i][2]);
}*/
/*
fixedPoints(ctx, "abc", ptmp);
drawCovector(ctx, cross(ptmp[0], ptmp[1]));
fixedPoints(ctx, "bca", ptmp);
drawCovector(ctx, cross(ptmp[0], ptmp[1]));
fixedPoints(ctx, "cab", ptmp);
drawCovector(ctx, cross(ptmp[0], ptmp[1]));
fixedPoints(ctx, "babcb", ptmp);
drawCovector(ctx, cross(ptmp[0], ptmp[1]));
fixedPoints(ctx, "cbcac", ptmp);
drawCovector(ctx, cross(ptmp[0], ptmp[1]));
fixedPoints(ctx, "acaba", ptmp);
drawCovector(ctx, cross(ptmp[0], ptmp[1]));
*/
cairo_restore(C);
}
void drawLimitCurve(DrawingContext *ctx)
{
int close_points = 0;
cairo_t *C = ctx->cairo;
cairo_save(C);
// cairo_set_source_rgb(C, 0.6, 0.6, 0.6);
cairo_set_source_rgb(C, 0, 0, 0);
if(ctx->limit_with_lines) {
int previous_inside = 0;
cairo_new_path(C);
for(int i = 0; i < ctx->limit_curve_count; i++) {
point_t p;
p.x = ctx->limit_curve[3*i];
p.y = ctx->limit_curve[3*i+1];
if(isInsideBB(ctx, p)) {
if(!previous_inside)
cairo_move_to(C, p.x, p.y);
else
cairo_line_to(C, p.x, p.y);
previous_inside = 1;
} else {
previous_inside = 0;
}
}
cairo_stroke(C);
} else {
for(int i = 0; i < ctx->limit_curve_count; i++) {
point_t p;
p.x = ctx->limit_curve[3*i];
p.y = ctx->limit_curve[3*i+1];
if(isInsideBB(ctx, p)) {
cairo_arc(C, p.x, p.y, 2.0/ctx->dim->scalefactor, 0, 2*M_PI);
cairo_close_path(C);
cairo_fill(C);
}
}
}
cairo_restore(C);
}
void drawCoxeterOrbit(DrawingContext *ctx)
{
gsl_matrix *rot = getTempMatrix(ctx->ws);
gsl_matrix **gen = getTempMatrices(ctx->ws, 6);
gsl_vector *eval = getTempVector(ctx->ws);
gsl_matrix *coxeter_fixedpoints = getTempMatrix(ctx->ws);
gsl_vector *startpoint_coxeterbasis = getTempVector(ctx->ws);
gsl_vector *startpoint_globalbasis = getTempVector(ctx->ws);
gsl_vector *startpoint_drawbasis = getTempVector(ctx->ws);
gsl_matrix **elements = getTempMatrices(ctx->ws, ctx->n_group_elements);
cairo_t *C = ctx->cairo;
vector_t cox[3][3];
vector_t abcb[3];
double ev[3];
vector_t v, start;
point_t p;
int first = 1;
cairo_save(C);
initializeTriangleGeneratorsCurrent(gen, ctx);
cairo_set_source_rgb(C, 0, 0, 1);
fixedPoints(ctx, "abc", cox[0]);
fixedPoints(ctx, "bca", cox[1]);
fixedPoints(ctx, "cab", cox[2]);
// fixedPoints(ctx, "babc", abcb);
wordEigenvalues(ctx, "abc", ev);
LOOP(i) LOOP(j) gsl_matrix_set(coxeter_fixedpoints, j, i, cox[0][i].x[j]);
initializeTriangleGeneratorsCurrent(gen, ctx);
gsl_matrix_set_identity(elements[0]);
for(int i = 1; i < ctx->n_group_elements; i++)
multiply(gen[ctx->group[i].letter], elements[ctx->group[i].parent->id], elements[i]);
//printf("coxeter eigenvalues: %f %f %f\n", ev[0], ev[1], ev[2]);
// LOOP(i) gsl_vector_set(startpoint_globalbasis, i, cox[1][0].x[i]);
gsl_vector_set(startpoint_drawbasis, 0, ctx->marking.x);
gsl_vector_set(startpoint_drawbasis, 1, ctx->marking.y);
gsl_vector_set(startpoint_drawbasis, 2, 1);
solve(ctx->cob, startpoint_drawbasis, startpoint_globalbasis, ctx->ws);
// solve(coxeter_fixedpoints, startpoint_globalbasis, startpoint_coxeterbasis, ctx->ws);
// LOOP(i) start.x[i] = gsl_vector_get(startpoint_coxeterbasis, i);
LOOP(i) start.x[i] = gsl_vector_get(startpoint_globalbasis, i);
/*
for(int t = -1000; t < 1000; t++) {
LOOP(i) v.x[i] = 0;
LOOP(i) LOOP(j) v.x[j] += pow(fabs(ev[i]),t*0.01)*start.x[i]*cox[0][i].x[j];
p = vectorToPoint(ctx, v);
if(first) {
cairo_move_to(C, p.x, p.y);
first = 0;
} else {
cairo_line_to(C, p.x, p.y);
}
}
cairo_stroke(C);
*/
for(int i = 0; i < ctx->n_group_elements; i++) {
v = apply(elements[i], start);
drawVector(ctx, v);
// ctx->limit_curve[3*i+2] = atan2(
// gsl_matrix_get(fixedpoints_pos, 2, column)/gsl_matrix_get(fixedpoints_pos, 0, column),
// gsl_matrix_get(fixedpoints_pos, 1, column)/gsl_matrix_get(fixedpoints_pos, 0, column));
}
/* for(int t = -20; t < 20; t++) {
LOOP(i) v.x[i] = 0;
LOOP(i) LOOP(j) v.x[j] += (ev[i]<0&&t%2?-1:1)*pow(fabs(ev[i]),t/3.0)*start.x[i]*cox[0][i].x[j];
drawVector(ctx, v);
}*/
cairo_set_source_rgb(C, 0, 0, 0);
// LOOP(i) drawVector(ctx, abcb[i]);
cairo_restore(C);
releaseTempMatrices(ctx->ws, 8 + ctx->n_group_elements);
releaseTempVectors(ctx->ws, 4);
}
void drawText(DrawingContext *ctx)
{
cairo_move_to(ctx->cairo, 15, 30);
cairo_set_source_rgb(ctx->cairo, 0, 0, 0);
char buf[100];
sprintf(buf, "t = exp(%.8f) = %.8f, s = exp(%.8f) = %.8f, marking = (%.5f, %.5f)", log(ctx->parameter), ctx->parameter, log(ctx->parameter2), ctx->parameter2, ctx->marking.x, ctx->marking.y);
cairo_show_text(ctx->cairo, buf);
}
// level 5: put everything together
void draw(DrawingContext *ctx)
{
cairo_t *C = ctx->cairo;
cairo_set_source_rgb(C, 1, 1, 1);
cairo_paint(C);
cairo_set_matrix(C, &ctx->dim->matrix);
// defaults; use save/restore whenever these are changed
cairo_set_line_width(C, 1.0/ctx->dim->scalefactor);
cairo_set_font_size(C, 16);
cairo_set_line_join(C, CAIRO_LINE_JOIN_BEVEL);
cairo_set_line_cap(C, CAIRO_LINE_CAP_ROUND);
if(ctx->limit_curve_count >= 0) {
if(ctx->show_limit)
drawLimitCurve(ctx);
if(ctx->show_dual_limit)
drawDualLimitCurve(ctx);
if(ctx->show_attractors)
drawAttractors(ctx);
if(ctx->show_rotated_reflectors)
drawRotatedReflectors(ctx);
if(ctx->show_reflectors)
drawReflectors(ctx);
if(ctx->show_boxes)
drawBoxes(ctx);
if(ctx->show_boxes2)
drawBoxes2(ctx);
if(ctx->show_marking)
{
cairo_set_source_rgb(C, 0, 0, 1);
drawPoint(ctx, ctx->marking);
}
if(ctx->show_coxeter_orbit)
drawCoxeterOrbit(ctx);
}
cairo_identity_matrix(C); // text is in screen coordinates
if(ctx->show_text)
drawText(ctx);
cairo_surface_flush(cairo_get_target(C));
}