naev 0.11.5
safelanes.c
Go to the documentation of this file.
1/*
2 * See Licensing and Copyright notice in naev.h
3 */
11#include <math.h>
12
13#if HAVE_SUITESPARSE_CHOLMOD_H
14#include <suitesparse/cholmod.h>
15#else /* HAVE_SUITESPARSE_CHOLMOD_H */
16#include <cholmod.h>
17#endif /* HAVE_SUITESPARSE_CHOLMOD_H */
18
19#if HAVE_OPENBLAS_CBLAS_H
20# include <openblas/cblas.h>
21#elif HAVE_CBLAS_OPENBLAS_H
22# include <cblas_openblas.h>
23#elif HAVE_CBLAS_HYPHEN_OPENBLAS_H
24# include <cblas-openblas.h>
25#elif HAVE_ACCELERATE_ACCELERATE_H
26# include <Accelerate/Accelerate.h>
27#elif HAVE_CBLAS_H
28# include <cblas.h>
29#elif HAVE_F77BLAS_H
30# include <f77blas.h>
31# define I_LOVE_FORTRAN 1
32#elif HAVE_OPENBLAS_F77BLAS_H
33# include <openblas/f77blas.h>
34# define I_LOVE_FORTRAN 1
35#endif
36
37#include "naev.h"
40#include "safelanes.h"
41
42#include "array.h"
43#include "conf.h"
44#include "log.h"
45#include "union_find.h"
46
47/*
48 * Global parameters.
49 */
50static const double ALPHA = 9.;
51static const double LAMBDA = 2e10;
52static const double JUMP_CONDUCTIVITY= 0.001;
53static const double MIN_ANGLE = M_PI/18.;
54enum {
58 SORTED = 1,
59 PACKED = 1,
61};
62
63/*
64 * Types.
65 */
69typedef enum VertexType_ {VERTEX_SPOB, VERTEX_JUMP} VertexType;
70
72typedef struct Vertex_ {
73 int system;
75 int index;
76} Vertex;
77
79typedef int Edge[2];
80
82typedef struct Faction_ {
83 int id;
84 double lane_length_per_presence;
85 double lane_base_cost;
86} Faction;
87
89typedef uint32_t FactionMask;
90static const FactionMask MASK_0 = 0, MASK_1 = 1;
91
92/*
93 * Global state.
94 */
95static cholmod_common C;
100static int *sys_to_first_edge;
102static int *lane_faction;
104static double **presence_budget;
105static int *tmp_spob_indices;
107static double *tmp_edge_conduct;
110static cholmod_triplet *stiff;
111static cholmod_sparse *QtQ;
112static cholmod_dense *ftilde;
113static cholmod_dense *utilde;
114static cholmod_dense **PPl;
115static double* cmp_key_ref;
118/*
119 * Prototypes.
120 */
121static int safelanes_buildOneTurn( int iters_done );
122static int safelanes_activateByGradient( cholmod_dense* Lambda_tilde, int iters_done );
123static void safelanes_initStacks (void);
124static void safelanes_initStacks_edge (void);
125static void safelanes_initStacks_faction (void);
126static void safelanes_initStacks_vertex (void);
127static void safelanes_initStacks_anchor (void);
128static void safelanes_initOptimizer (void);
129static void safelanes_destroyOptimizer (void);
130static void safelanes_destroyStacks (void);
131static void safelanes_destroyTmp (void);
132static void safelanes_initStiff (void);
133static double safelanes_initialConductivity ( int ei );
134static void safelanes_updateConductivity ( int ei_activated );
135static void safelanes_initQtQ (void);
136static void safelanes_initFTilde (void);
137static void safelanes_initPPl (void);
138static int safelanes_triangleTooFlat( const vec2* m, const vec2* n, const vec2* p, double lmn );
139static int vertex_faction( int vi );
140static const vec2* vertex_pos( int vi );
141static inline int FACTION_ID_TO_INDEX( int id );
142static inline FactionMask MASK_ANY_FACTION();
143static inline FactionMask MASK_ONE_FACTION( int id );
144static inline FactionMask MASK_COMPROMISE( int id1, int id2 );
145static int cmp_key( const void* p1, const void* p2 );
146static inline void triplet_entry( cholmod_triplet* m, int i, int j, double v );
147static cholmod_dense* safelanes_sliceByPresence( cholmod_dense* m, double* sysPresence );
148static cholmod_dense* ncholmod_ddmult( cholmod_dense* A, int transA, cholmod_dense* B );
149static double safelanes_row_dot_row( cholmod_dense* A, cholmod_dense* B, int i, int j );
150
154static inline void array_push_back_edge( Edge **a, int v0, int v1 )
155{
156 int *vs = array_grow( a );
157 vs[0] = v0;
158 vs[1] = v1;
159}
160
164void safelanes_init (void)
165{
166 cholmod_start( &C );
167 /* Ideally we would want to recalculate here, but since we load the first
168 * save and try to use unidiffs there, we instead defer the safe lane
169 * computation to only if necessary after loading save unidiffs. */
170 /* safelanes_recalculate(); */
171}
172
177{
180 cholmod_finish( &C );
181}
182
190SafeLane* safelanes_get( int faction, int standing, const StarSystem* system )
191{
193
194 for (int i=sys_to_first_edge[system->id]; i<sys_to_first_edge[1+system->id]; i++) {
195 SafeLane *l;
196 const Vertex *v[2];
197 int lf = lane_faction[i];
198
199 /* No lane on edge. */
200 if (lf <= 0)
201 continue;
202
203 /* Filter by standing. */
204 if (faction >= 0) {
205 if (standing==0) {
206 /* Only match exact faction. */
207 if (lf != faction)
208 continue;
209 }
210 else {
211 /* Try to do more advanced matching. */
212 int fe = areEnemies(faction,lf);
213 int fa = areAllies(faction,lf);
214 int skip = 1;
215 if ((standing & SAFELANES_FRIENDLY) && fa)
216 skip = 0;
217 if ((standing & SAFELANES_NEUTRAL) && !fe)
218 skip = 0;
219 if ((standing & SAFELANES_HOSTILE) && fe)
220 skip = 0;
221 if (skip)
222 continue;
223 }
224 }
225
226 for (int j=0; j<2; j++)
227 v[j] = &vertex_stack[edge_stack[i][j]];
228
229 l = &array_grow( &out );
230 l->faction = lane_faction[i];
231 for (int j=0; j<2; j++) {
232 switch (v[j]->type) {
233 case VERTEX_SPOB:
234 l->point_type[j] = SAFELANE_LOC_SPOB;
235 l->point_id[j] = system->spobs[v[j]->index]->id;
236 break;
237 case VERTEX_JUMP:
238 l->point_type[j] = SAFELANE_LOC_DEST_SYS;
239 l->point_id[j] = system->jumps[v[j]->index].targetid;
240 break;
241 default:
242 ERR( _("Safe-lane vertex type is invalid.") );
243 }
244 }
245 }
246 return out;
247}
248
253{
254 Uint32 time = SDL_GetTicks();
255
256 /* Don't recompute on exit. */
257 if (naev_isQuit())
258 return;
259
262 for (int iters_done=0; safelanes_buildOneTurn(iters_done) > 0; iters_done++)
263 ;
265 /* Stacks remain available for queries. */
266 time = SDL_GetTicks() - time;
267 if (conf.devmode)
268 DEBUG( n_("Charted safe lanes for %d object in %.3f s", "Charted safe lanes for %d objects in %.3f s", array_size(vertex_stack)), array_size(vertex_stack), time/1000. );
269
271}
272
277{
279}
280
292
297{
298 for (int i=0; i<array_size(PPl); i++)
299 cholmod_free_dense( &PPl[i], &C );
300 array_free( PPl );
301 PPl = NULL;
302 cholmod_free_dense( &utilde, &C ); /* CAUTION: if we instead save it, ensure it's updated after the final activateByGradient. */
303 cholmod_free_dense( &ftilde, &C );
304 cholmod_free_sparse( &QtQ, &C );
305 cholmod_free_triplet( &stiff, &C );
306}
307
311static int safelanes_buildOneTurn( int iters_done )
312{
313 cholmod_sparse *stiff_s;
314 cholmod_factor *stiff_f;
315 cholmod_dense *_QtQutilde, *Lambda_tilde, *Y_workspace, *E_workspace;
316 int turns_next_time;
317 double zero[] = {0, 0}, neg_1[] = {-1, 0};
318
319 Y_workspace = E_workspace = Lambda_tilde = NULL;
320 stiff_s = cholmod_triplet_to_sparse( stiff, 0, &C );
321 stiff_f = cholmod_analyze( stiff_s, &C );
322 cholmod_factorize( stiff_s, stiff_f, &C );
323 cholmod_solve2( CHOLMOD_A, stiff_f, ftilde, NULL, &utilde, NULL, &Y_workspace, &E_workspace, &C );
324 _QtQutilde = cholmod_zeros( utilde->nrow, utilde->ncol, CHOLMOD_REAL, &C );
325 cholmod_sdmult( QtQ, 0, neg_1, zero, utilde, _QtQutilde, &C );
326 cholmod_solve2( CHOLMOD_A, stiff_f, _QtQutilde, NULL, &Lambda_tilde, NULL, &Y_workspace, &E_workspace, &C );
327 cholmod_free_dense( &_QtQutilde, &C );
328 cholmod_free_dense( &Y_workspace, &C );
329 cholmod_free_dense( &E_workspace, &C );
330 cholmod_free_factor( &stiff_f, &C );
331 cholmod_free_sparse( &stiff_s, &C );
332 turns_next_time = safelanes_activateByGradient( Lambda_tilde, iters_done );
333 cholmod_free_dense( &Lambda_tilde, &C );
334
335 return turns_next_time;
336}
337
341static void safelanes_initStacks (void)
342{
344 safelanes_initStacks_faction(); /* Dependency for vertex. */
345 safelanes_initStacks_vertex(); /* Dependency for edge. */
348}
349
354{
355 const StarSystem *systems_stack = system_getAll();
356
362 for (int system=0; system<array_size(systems_stack); system++) {
363 const StarSystem *sys = &systems_stack[system];
364 if (sys_isFlag( sys, SYSTEM_NOLANES )) {
366 continue;
367 }
368
369 for (int i=0; i<array_size(sys->spobs); i++) {
370 const Spob *p = sys->spobs[i];
371 if (spob_isFlag( p, SPOB_NOLANES ))
372 continue;
373 if (p->presence.base!=0. || p->presence.bonus!=0.) {
374 Vertex v = {.system = system, .type = VERTEX_SPOB, .index = i};
377 }
378 }
379
380 for (int i=0; i<array_size(sys->jumps); i++) {
381 const JumpPoint *jp = &sys->jumps[i];
382 if (jp_isFlag( jp, JP_HIDDEN | JP_EXITONLY | JP_NOLANES ))
383 continue;
384 Vertex v = {.system = system, .type = VERTEX_JUMP, .index = i};
386 if (jp->targetid < system && jp->returnJump != NULL)
387 for (int j=sys_to_first_vertex[jp->targetid]; j < sys_to_first_vertex[1+jp->targetid]; j++)
388 if (vertex_stack[j].type == VERTEX_JUMP && jp->returnJump == &jp->target->jumps[vertex_stack[j].index]) {
390 break;
391 }
392 }
393
395 }
396 //array_shrink( &vertex_stack );
397 //array_shrink( &sys_to_first_vertex );
398
399 vertex_fmask = calloc( array_size(vertex_stack), sizeof(FactionMask) );
400}
401
406{
407 const StarSystem *systems_stack = system_getAll();
408
413 tmp_edge_conduct = array_create( double );
414 for (int system=0; system<array_size(systems_stack); system++) {
415 for (int i=sys_to_first_vertex[system]; i < sys_to_first_vertex[1+system]; i++) {
416 const vec2 *pi = vertex_pos( i );
417 for (int j=sys_to_first_vertex[system]; j < i; j++) {
418 const vec2 *pj = vertex_pos( j );
419 double lij = vec2_dist( pi, pj );
420 int has_approx_midpoint = 0;
421 for (int k=sys_to_first_vertex[system]; k < sys_to_first_vertex[1+system]; k++)
422 if (k!=i && k!=j && safelanes_triangleTooFlat( pi, pj, vertex_pos( k ), lij )) {
423 has_approx_midpoint = 1;
424 break;
425 }
426 if (has_approx_midpoint)
427 continue; /* The edge from i to j is disallowed in favor of a path through k. */
431 }
432 }
434 }
437
440 memset( lane_faction, 0, array_size(lane_faction)*sizeof(lane_faction[0]) );
441}
442
447{
448 int *faction_all;
449 const StarSystem *systems_stack;
450
452 faction_all = faction_getAllVisible();
453 for (int fi=0; fi<array_size(faction_all); fi++) {
454 int f = faction_all[fi];
455 Faction rec = {.id = f, .lane_length_per_presence = faction_lane_length_per_presence(f), .lane_base_cost = faction_lane_base_cost(f)};
456 if (rec.lane_length_per_presence > 0.)
458 }
459 array_free( faction_all );
461 assert( "FactionMask size is sufficient" && (size_t)array_size(faction_stack) <= 8*sizeof(FactionMask) );
462
465 for (int fi=0; fi<array_size(faction_stack); fi++) {
467 for (int s=0; s<array_size(systems_stack); s++) {
468 const StarSystem *sys = &systems_stack[s];
469 double budget = system_getPresence( sys, faction_stack[fi].id );
470 array_push_back( &presence_budget[fi], budget );
471 }
472 }
473}
474
481{
482 int *anchor_systems;
483 int nsys = array_size(sys_to_first_vertex) - 1;
484 unionfind_init( &tmp_sys_uf, nsys );
485 for (int i=0; i<array_size(tmp_jump_edges); i++)
487 anchor_systems = unionfind_findall( &tmp_sys_uf );
488 tmp_anchor_vertices = array_create_size( int, array_size(anchor_systems) );
489
490 /* Add an anchor vertex per system, but only if there actually is a vertex in the system. */
491 for (int i=0; i<array_size(anchor_systems); i++)
492 if (sys_to_first_vertex[anchor_systems[i]] < sys_to_first_vertex[1+anchor_systems[i]])
494 array_free( anchor_systems );
495}
496
500static void safelanes_destroyStacks (void)
501{
504 vertex_stack = NULL;
505 free( vertex_fmask );
506 vertex_fmask = NULL;
508 sys_to_first_vertex = NULL;
510 edge_stack = NULL;
512 sys_to_first_edge = NULL;
513 for (int i=0; i<array_size(presence_budget); i++)
516 presence_budget = NULL;
518 faction_stack = NULL;
520 lane_faction = NULL;
522 lane_fmask = NULL;
523}
524
540
544static void safelanes_initStiff (void)
545{
546 int nnz, v;
547 double max_conductivity;
548
549 cholmod_free_triplet( &stiff, &C );
552 stiff = cholmod_allocate_triplet( v, v, nnz, STORAGE_MODE_UPPER_TRIANGULAR_PART, CHOLMOD_REAL, &C );
553 /* Populate triplets: internal edges (ii ij jj), implicit jump connections (ditto), anchor conditions. */
554 for (int i=0; i<array_size(edge_stack); i++) {
555 triplet_entry( stiff, edge_stack[i][0], edge_stack[i][0], +tmp_edge_conduct[i] );
556 triplet_entry( stiff, edge_stack[i][0], edge_stack[i][1], -tmp_edge_conduct[i] );
557 triplet_entry( stiff, edge_stack[i][1], edge_stack[i][1], +tmp_edge_conduct[i] );
558 }
559 for (int i=0; i<array_size(tmp_jump_edges); i++) {
560 triplet_entry( stiff, tmp_jump_edges[i][0], tmp_jump_edges[i][0], +JUMP_CONDUCTIVITY );
561 triplet_entry( stiff, tmp_jump_edges[i][0], tmp_jump_edges[i][1], -JUMP_CONDUCTIVITY );
562 triplet_entry( stiff, tmp_jump_edges[i][1], tmp_jump_edges[i][1], +JUMP_CONDUCTIVITY );
563 }
564 /* Add a Robin boundary condition, using the max conductivity (after activation) for spectral reasons. */
565 max_conductivity = JUMP_CONDUCTIVITY/(1+ALPHA);
566 for (int i=0; i<array_size(edge_stack); i++)
567 max_conductivity = MAX( max_conductivity, tmp_edge_conduct[i] );
568 max_conductivity = MAX( JUMP_CONDUCTIVITY, (1+ALPHA)*max_conductivity ); /* Activation scales entries by 1+ALPHA later. */
569 for (int i=0; i<array_size(tmp_anchor_vertices); i++)
570 triplet_entry( stiff, tmp_anchor_vertices[i], tmp_anchor_vertices[i], max_conductivity );
571#if DEBUGGING
572 assert( stiff->nnz == stiff->nzmax );
573 assert( cholmod_check_triplet( stiff, &C) );
574#endif /* DEBUGGING */
575}
576
582static double safelanes_initialConductivity ( int ei )
583{
584 double *sv = stiff->x;
585 return lane_faction[ei] ? sv[3*ei]/(1+ALPHA) : sv[3*ei];
586}
587
592static void safelanes_updateConductivity ( int ei_activated )
593{
594 double *sv = stiff->x;
595 for (int i=3*ei_activated; i<3*(ei_activated+1); i++)
596 sv[i] *= 1+ALPHA;
597}
598
602static void safelanes_initQtQ (void)
603{
604 cholmod_sparse *Q;
605
606 cholmod_free_sparse( &QtQ, &C );
607 /* Form Q, the edge-vertex projection where (Dirac notation) Q |edge> = |edge[0]> - |edge[1]>. It has a +1 and -1 per column. */
608 Q = cholmod_allocate_sparse( array_size(vertex_stack), array_size(edge_stack), 2*array_size(edge_stack),
609 SORTED, PACKED, STORAGE_MODE_UNSYMMETRIC, CHOLMOD_REAL, &C );
610 ((int*)Q->p)[0] = 0;
611 for (int i=0; i<array_size(edge_stack); i++) {
612 ((int*)Q->p)[i+1] = 2*(i+1);
613 ((int*)Q->i)[2*i+0] = edge_stack[i][0];
614 ((int*)Q->i)[2*i+1] = edge_stack[i][1];
615 ((double*)Q->x)[2*i+0] = +1;
616 ((double*)Q->x)[2*i+1] = -1;
617 }
618#if DEBUGGING
619 assert( cholmod_check_sparse( Q, &C ) );
620#endif /* DEBUGGING */
621 QtQ = cholmod_aat( Q, NULL, 0, MODE_NUMERICAL, &C );
622 cholmod_free_sparse( &Q, &C );
623}
624
628static void safelanes_initFTilde (void)
629{
630 cholmod_sparse *eye = cholmod_speye( array_size(vertex_stack), array_size(vertex_stack), CHOLMOD_REAL, &C );
631 cholmod_sparse *sp = cholmod_submatrix( eye, NULL, -1, tmp_spob_indices, array_size(tmp_spob_indices), 1, SORTED, &C );
632 cholmod_free_dense( &ftilde, &C );
633 ftilde = cholmod_sparse_to_dense( sp, &C );
634 cholmod_free_sparse( &sp, &C );
635 cholmod_free_sparse( &eye, &C );
636}
637
641static void safelanes_initPPl (void)
642{
643 int *component;
645
646 for (int fi=0; fi<array_size(PPl); fi++)
647 cholmod_free_dense( &PPl[fi], &C );
648 array_free( PPl );
649 PPl = array_create_size( cholmod_dense*, array_size(faction_stack) );
650 for (int fi=0; fi<array_size(faction_stack); fi++)
651 array_push_back( &PPl, cholmod_zeros( np, np, CHOLMOD_REAL, &C ) );
652
653 /* Form P, the pair-vertex projection where (Dirac notation) P |pair(i,j)> = |i> - |j>. It has a +1 and -1 per column. */
654 /* At least, pretend we did. We want (PD)(PD)*, where D is a diagonal matrix whose pair(i,j) are these presence sums: */
655
656 component = calloc( np, sizeof(int) );
657 for (int i=0; i<np; i++)
658 component[i] = unionfind_find( &tmp_sys_uf, vertex_stack[tmp_spob_indices[i]].system );
659
660 for (int i=0; i<np; i++) {
661 double *Di;
663 Spob *pnt = system_getIndex( sys )->spobs[vertex_stack[tmp_spob_indices[i]].index];
664 double pres = pnt->presence.base + pnt->presence.bonus; /* TODO distinguish between base and bonus? */
665 int fi = FACTION_ID_TO_INDEX( pnt->presence.faction );
666 if (fi < 0)
667 continue;
668 Di = PPl[fi]->x;
669 for (int j=0; j<i; j++)
670 if (component[i] == component[j])
671 Di[np*i+j] += pres;
672 for (int j=i+1; j<np; j++)
673 if (component[i] == component[j])
674 Di[np*j+i] += pres;
675 }
676
677 /* At this point, PPl[fi]->x[np*i+j] holds the pair(i,j) entry of D. */
678 for (int fi=0; fi<array_size(faction_stack); fi++)
679 for (int i=0; i<np; i++)
680 for (int j=0; j<i; j++) {
681 double d = ((double*)PPl[fi]->x)[np*i+j];
682 d *= d;
683 ((double*)PPl[fi]->x)[np*i+j] = -d;
684 ((double*)PPl[fi]->x)[np*j+i] = -d;
685 ((double*)PPl[fi]->x)[np*i+i] += d;
686 ((double*)PPl[fi]->x)[np*j+j] += d;
687 }
688
689 free( component );
690}
691
696static int safelanes_activateByGradient( cholmod_dense* Lambda_tilde, int iters_done )
697{
698 int *facind_opts, *edgeind_opts, turns_next_time;
699 double *facind_vals, Linv;
700 cholmod_dense **lal;
701 size_t *lal_bases, lal_base;
703 lal = calloc( array_size(faction_stack), sizeof(cholmod_dense*) );
704 lal_bases = calloc( array_size(faction_stack), sizeof(size_t) );
705 edgeind_opts = array_create( int );
706 facind_opts = array_create_size( int, array_size(faction_stack) );
707 facind_vals = array_create_size( double, array_size(faction_stack) );
708 for (int fi=0; fi<array_size(faction_stack); fi++) {
709 array_push_back( &facind_opts, fi );
710 array_push_back( &facind_vals, 0 );
711 }
712 turns_next_time = 0;
713
714 for (int si=0; si<array_size(sys_to_first_vertex)-1; si++) {
715 /* Factions with most presence here choose first. */
716 StarSystem *sys = system_getIndex( si );
717 for (int fi=0; fi<array_size(faction_stack); fi++)
718 facind_vals[fi] = -system_getPresence( sys, faction_stack[fi].id ); /* FIXME: Is this better, or presence_budget? */
719 cmp_key_ref = facind_vals;
720 qsort( facind_opts, array_size(faction_stack), sizeof(int), cmp_key );
721
722 for (int fii=0; fii<array_size(faction_stack); fii++) {
723 int ei_best;
724 double cost_best, cost_cheapest_other;
725 int fi = facind_opts[fii];
726 size_t sys_base = sys_to_first_vertex[si];
727 double score_best = 0.; /* Negative scores get ignored. */
728
729 /* Get the base index to use for this system. Save the value we expect to be the next iteration's base index.
730 * The current system's rows are in the lal[fi] matrix if there's presence at the time we slice it.
731 * What we know is whether there's presence *now*. We can use that as a proxy and fix lal_bases[fi] if we
732 * deplete this presence before constructing lal[fi]. This is tricky, so there are assertions below,
733 * which can warn us if we fuck this up. */
734 lal_base = lal_bases[fi];
735 if (presence_budget[fi][si] <= 0.)
736 continue;
737 /* We "should" find these DoF's interesting if/when we slice, and will unless we deplete this presence first. */
738 lal_bases[fi] += sys_to_first_vertex[1+si] - sys_to_first_vertex[si];
739
740 array_resize( &edgeind_opts, 0 );
741 for (int ei=sys_to_first_edge[si]; ei<sys_to_first_edge[1+si]; ei++) {
742 int sis = edge_stack[ei][0];
743 int sjs = edge_stack[ei][1];
744 int disconnecting = iters_done && !(vertex_fmask[sis] & (MASK_1<<fi)) && !(vertex_fmask[sjs] & (MASK_1<<fi));
746 if (!lane_faction[ei]
747 && !disconnecting
748 && presence_budget[fi][si] >= cost
749 && (lane_fmask[ei] & (MASK_1<<fi)))
750 array_push_back( &edgeind_opts, ei );
751 }
752
753 if (array_size(edgeind_opts) == 0) {
754 presence_budget[fi][si] = 0.; /* Nothing to build here! Tell ourselves to stop trying. */
755 if (lal[fi] == NULL)
756 lal_bases[fi] -= sys_to_first_vertex[1+si] - sys_to_first_vertex[si];
757 continue;
758 }
759
760 ei_best = edgeind_opts[0];
762 cost_cheapest_other = +HUGE_VAL;
763 if (array_size(edgeind_opts) > 0) {
764 /* There's an actual choice. Search for the best option. Lower is better. */
765 for (int eii=0; eii<array_size(edgeind_opts); eii++) {
766 int ei = edgeind_opts[eii];
767 int sis = edge_stack[ei][0];
768 int sjs = edge_stack[ei][1];
769 double score = 0.;
770 double cost;
771
772 if (lal[fi] == NULL) { /* Is it time to evaluate the lazily-calculated matrix? */
773 cholmod_dense *lamt = safelanes_sliceByPresence( Lambda_tilde, presence_budget[fi] );
774 lal[fi] = ncholmod_ddmult( lamt, 0, PPl[fi] );
775 cholmod_free_dense( &lamt, &C );
776 }
777
778 /* Evaluate (LUTll[0,0] + LUTll[1,1] - LUTll[0,1] - LUTll[1,0]), */
779 /* where LUTll = np.dot( lal[[sis,sjs],:] , utilde[[sis,sjs],:].T ) */
780 score += safelanes_row_dot_row( utilde, lal[fi], sis, sis - sys_base + lal_base );
781 score += safelanes_row_dot_row( utilde, lal[fi], sjs, sjs - sys_base + lal_base );
782 score -= safelanes_row_dot_row( utilde, lal[fi], sjs, sis - sys_base + lal_base );
783 score -= safelanes_row_dot_row( utilde, lal[fi], sis, sjs - sys_base + lal_base );
785 score *= ALPHA * Linv * Linv;
786 score += LAMBDA;
787
789 if (score < score_best) {
790 ei_best = ei;
791 score_best = score;
792 cost_cheapest_other = MIN( cost_cheapest_other, cost_best );
793 cost_best = cost;
794 }
795 else
796 cost_cheapest_other = MIN( cost_cheapest_other, cost );
797 }
798 }
799
800 /* Ignore positive scores. */
801 if (score_best >= 0.)
802 continue;
803
804 /* Add the lane. */
805 presence_budget[fi][si] -= cost_best;
806 if (presence_budget[fi][si] >= cost_cheapest_other)
807 turns_next_time++;
808 else {
809 presence_budget[fi][si] = 0.; /* Nothing more to do here; tell ourselves. */
810 if (lal[fi] == NULL)
811 lal_bases[fi] -= sys_to_first_vertex[1+si] - sys_to_first_vertex[si];
812 }
814 vertex_fmask[edge_stack[ei_best][0]] |= (MASK_1<<fi);
815 vertex_fmask[edge_stack[ei_best][1]] |= (MASK_1<<fi);
816 lane_faction[ ei_best ] = faction_stack[fi].id;
817 }
818 }
819
820#if DEBUGGING
821 for (int fi=0; fi<array_size(faction_stack); fi++)
822 if (lal[fi] != NULL)
823 assert( "Correctly tracked row offsets between the 'lal' and 'utilde' matrices" && lal[fi]->nrow == lal_bases[fi] );
824#endif /* DEBUGGING */
825
826 for (int fi=0; fi<array_size(faction_stack); fi++)
827 cholmod_free_dense( &lal[fi], &C );
828 free( lal );
829 free( lal_bases );
830 array_free( edgeind_opts );
831 array_free( facind_vals );
832 array_free( facind_opts );
833
834 return turns_next_time;
835}
836
838static int cmp_key( const void* p1, const void* p2 )
839{
840 double d = cmp_key_ref[*(int*)p1] - cmp_key_ref[*(int*)p2];
841 return SIGN( d );
842}
843
847static int safelanes_triangleTooFlat( const vec2* m, const vec2* n, const vec2* p, double lmn )
848{
849 const double MAX_COSINE = cos(MIN_ANGLE);
850 double lnp = vec2_dist( n, p );
851 double lmp = vec2_dist( m, p );
852 double dpn = ((n->x-m->x)*(n->x-p->x) + (n->y-m->y)*(n->y-p->y)) / ( lmn * lnp );
853 double dpm = ((m->x-n->x)*(m->x-p->x) + (m->y-n->y)*(m->y-p->y)) / ( lmn * lmp );
854 return (dpn > MAX_COSINE && lnp < lmn) || (dpm > MAX_COSINE && lmp < lmn);
855}
856
860static int vertex_faction( int vi )
861{
862 const StarSystem *sys = system_getIndex(vertex_stack[vi].system);
863 switch (vertex_stack[vi].type) {
864 case VERTEX_SPOB:
865 return sys->spobs[vertex_stack[vi].index]->presence.faction;
866 case VERTEX_JUMP:
867 return -1;
868 default:
869 ERR( _("Safe-lane vertex type is invalid.") );
870 }
871}
872
876static const vec2* vertex_pos( int vi )
877{
878 const StarSystem *sys = system_getIndex(vertex_stack[vi].system);
879 switch (vertex_stack[vi].type) {
880 case VERTEX_SPOB:
881 return &sys->spobs[vertex_stack[vi].index]->pos;
882 case VERTEX_JUMP:
883 return &sys->jumps[vertex_stack[vi].index].pos;
884 default:
885 ERR( _("Safe-lane vertex type is invalid.") );
886 }
887}
888
890static inline int FACTION_ID_TO_INDEX( int id )
891{
892 for (int i=0; i<array_size(faction_stack) && faction_stack[i].id <= id; i++)
893 if (faction_stack[i].id == id)
894 return i;
895 return -1;
896}
897
900{
901 return ~MASK_0;
902}
903
905static inline FactionMask MASK_ONE_FACTION( int id )
906{
907 int ind = FACTION_ID_TO_INDEX( id );
908 return ind>0 ? (MASK_1)<<ind : MASK_ANY_FACTION();
909}
910
912static inline FactionMask MASK_COMPROMISE( int id1, int id2 )
913{
914 FactionMask m1 = MASK_ONE_FACTION(id1), m2 = MASK_ONE_FACTION(id2);
915 return (m1 & m2) ? (m1 & m2) : (m1 | m2); /* Any/Any -> any, Any/f -> just f, f1/f2 -> either. */
916}
917
918static inline void triplet_entry( cholmod_triplet* m, int i, int j, double x )
919{
920 ((int*)m->i)[m->nnz] = i;
921 ((int*)m->j)[m->nnz] = j;
922 ((double*)m->x)[m->nnz] = x;
923 m->nnz++;
924}
925
929static cholmod_dense* safelanes_sliceByPresence( cholmod_dense* m, double* sysPresence )
930{
931 size_t nr, nc, in_r, out_r;
932 cholmod_dense *out;
933
934 nr = 0;
935 nc = m->ncol;
936 for (int si=0; si < array_size(sys_to_first_vertex) - 1; si++)
937 if (sysPresence[si] > 0)
939
940 out = cholmod_allocate_dense( nr, nc, nr, CHOLMOD_REAL, &C );
941
942 in_r = out_r = 0;
943 for (int si=0; si < array_size(sys_to_first_vertex) - 1; si++) {
944 int sz = sys_to_first_vertex[1+si] - sys_to_first_vertex[si];
945 if (sysPresence[si] > 0) {
946 for (size_t c = 0; c < nc; c++)
947 memcpy( &((double*)out->x)[c*out->d + out_r], &((double*)m->x)[c*m->d + in_r], sz * sizeof(double) );
948 out_r += sz;
949 }
950 in_r += sz;
951 }
952 return out;
953}
954
956static cholmod_dense* ncholmod_ddmult( cholmod_dense* A, int transA, cholmod_dense* B )
957{
958#if I_LOVE_FORTRAN
959 blasint M = transA ? A->ncol : A->nrow, K = transA ? A->nrow : A->ncol, N = B->ncol, lda = A->d, ldb = B->d, ldc = M;
960 assert( K == (blasint) B->nrow );
961 cholmod_dense *out = cholmod_allocate_dense( M, N, M, CHOLMOD_REAL, &C );
962 double alpha = 1., beta = 0.;
963 BLASFUNC(dgemm)( transA ? "T" : "N", "N", &M, &N, &K, &alpha, A->x, &lda, B->x, &ldb, &beta, out->x, &ldc);
964#else /* I_LOVE_FORTRAN */
965 size_t M = transA ? A->ncol : A->nrow, K = transA ? A->nrow : A->ncol, N = B->ncol;
966 assert( K == B->nrow );
967 cholmod_dense *out = cholmod_allocate_dense( M, N, M, CHOLMOD_REAL, &C );
968 cblas_dgemm( CblasColMajor, transA?CblasTrans:CblasNoTrans, CblasNoTrans, M, N, K, 1, A->x, A->d, B->x, B->d, 0, out->x, out->d);
969#endif /* I_LOVE_FORTRAN */
970 return out;
971}
972
974static double safelanes_row_dot_row( cholmod_dense* A, cholmod_dense* B, int i, int j )
975{
976 assert( A->ncol == B->ncol );
977#if I_LOVE_FORTRAN
978 blasint N = A->ncol, incA = A->d, incB = B->d;
979 return BLASFUNC(ddot)( &N, (double*)A->x + i, &incA, (double*)B->x + j, &incB);
980#else /* I_LOVE_FORTRAN */
981 return cblas_ddot( A->ncol, (double*)A->x + i, A->d, (double*)B->x + j, B->d);
982#endif /* I_LOVE_FORTRAN */
983}
Provides macros to work with dynamic arrays.
#define array_free(ptr_array)
Frees memory allocated and sets array to NULL.
Definition array.h:158
#define array_resize(ptr_array, new_size)
Resizes the array to accomodate new_size elements.
Definition array.h:112
#define array_create_size(basic_type, capacity)
Creates a new dynamic array of ‘basic_type’ with an initial capacity.
Definition array.h:102
static ALWAYS_INLINE int array_size(const void *array)
Returns number of elements in the array.
Definition array.h:168
#define array_grow(ptr_array)
Increases the number of elements by one and returns the last element.
Definition array.h:119
#define array_shrink(ptr_array)
Shrinks memory to fit only ‘size’ elements.
Definition array.h:149
#define array_push_back(ptr_array, element)
Adds a new element at the end of the array.
Definition array.h:129
#define array_create(basic_type)
Creates a new dynamic array of ‘basic_type’.
Definition array.h:93
StarSystem * systems_stack
Definition space.c:92
int areEnemies(int a, int b)
Checks whether two factions are enemies.
Definition faction.c:1227
double faction_lane_base_cost(int f)
Gets the faction's weight for patrolled safe-lane construction;.
Definition faction.c:438
double faction_lane_length_per_presence(int f)
Gets the faction's weight for patrolled safe-lane construction (0 means they don't build lanes).
Definition faction.c:426
int * faction_getAllVisible(void)
Returns all non-invisible faction IDs in an array (array.h).
Definition faction.c:207
int areAllies(int a, int b)
Checks whether two factions are allies or not.
Definition faction.c:1253
int naev_isQuit(void)
Get if Naev is trying to quit.
Definition naev.c:158
Header file with generic functions and naev-specifics.
#define MIN(x, y)
Definition naev.h:40
#define SIGN(x)
Definition naev.h:43
#define MAX(x, y)
Definition naev.h:39
static const double c[]
Definition rng.c:264
static const double d[]
Definition rng.c:273
static int safelanes_calculated_once
Definition safelanes.c:116
static FactionMask MASK_COMPROMISE(int id1, int id2)
A mask with appropriate lane-building rights given one faction ID owning each endpoint.
Definition safelanes.c:912
@ STORAGE_MODE_UPPER_TRIANGULAR_PART
Definition safelanes.c:57
@ PACKED
Definition safelanes.c:59
@ SORTED
Definition safelanes.c:58
@ MODE_NUMERICAL
Definition safelanes.c:60
@ STORAGE_MODE_LOWER_TRIANGULAR_PART
Definition safelanes.c:55
@ STORAGE_MODE_UNSYMMETRIC
Definition safelanes.c:56
static FactionMask * vertex_fmask
Definition safelanes.c:97
static cholmod_dense * safelanes_sliceByPresence(cholmod_dense *m, double *sysPresence)
Construct the matrix-slice of m, selecting those rows where the corresponding presence value is posit...
Definition safelanes.c:929
static void array_push_back_edge(Edge **a, int v0, int v1)
Like array_push_back( a, Edge{v0, v1} ), but achievable in C. :-P.
Definition safelanes.c:154
static FactionMask * lane_fmask
Definition safelanes.c:103
static const double MIN_ANGLE
Definition safelanes.c:53
static void safelanes_initStiff(void)
Sets up the stiffness matrix.
Definition safelanes.c:544
static void safelanes_destroyTmp(void)
Tears down the local faction/object stacks.
Definition safelanes.c:528
static const double ALPHA
Definition safelanes.c:50
static cholmod_dense * ftilde
Definition safelanes.c:112
static void safelanes_updateConductivity(int ei_activated)
Updates the stiffness matrix to account for the given edge being activated.
Definition safelanes.c:592
static void safelanes_destroyStacks(void)
Tears down the local faction/object stacks.
Definition safelanes.c:500
int safelanes_calculated(void)
Whether or not the safe lanes have been calculated at least once.
Definition safelanes.c:276
void safelanes_destroy(void)
Shuts down the safelanes system.
Definition safelanes.c:176
static double ** presence_budget
Definition safelanes.c:104
static Faction * faction_stack
Definition safelanes.c:101
static Edge * edge_stack
Definition safelanes.c:99
static void safelanes_initStacks(void)
Sets up the local faction/object stacks.
Definition safelanes.c:341
static int safelanes_activateByGradient(cholmod_dense *Lambda_tilde, int iters_done)
Per-system, per-faction, activates the affordable lane with best (grad phi)/L.
Definition safelanes.c:696
static void safelanes_initStacks_anchor(void)
Identifies anchor points: The universe graph (with in-system and 2-way-jump edges) could have many co...
Definition safelanes.c:480
static void safelanes_initPPl(void)
Sets up the PPl matrices appearing in the gradient formula.
Definition safelanes.c:641
void safelanes_init(void)
Initializes the safelanes system.
Definition safelanes.c:164
static int vertex_faction(int vi)
Return the vertex's owning faction (ID, not faction_stack index), or -1 if not applicable.
Definition safelanes.c:860
static const vec2 * vertex_pos(int vi)
Return the vertex's coordinates within its system (by reference since our vec2's are fat).
Definition safelanes.c:876
static FactionMask MASK_ANY_FACTION()
Return a mask matching any faction.
Definition safelanes.c:899
static double safelanes_initialConductivity(int ei)
Returns the initial conductivity value (1/length) for edge ei. The live value is stored in the stiffn...
Definition safelanes.c:582
static cholmod_dense ** PPl
Definition safelanes.c:114
static void safelanes_initStacks_edge(void)
Sets up the local stacks with entry per edge. Faction stack must be set up.
Definition safelanes.c:405
static UnionFind tmp_sys_uf
Definition safelanes.c:109
static double safelanes_row_dot_row(cholmod_dense *A, cholmod_dense *B, int i, int j)
Return the i,j entry of A*B', or equivalently the dot product of row i of A with row j of B.
Definition safelanes.c:974
static int cmp_key(const void *p1, const void *p2)
It's a qsort comparator. Set the cmp_key_ref pointer prior to use, or else.
Definition safelanes.c:838
static cholmod_triplet * stiff
Definition safelanes.c:110
static int safelanes_buildOneTurn(int iters_done)
Run a round of optimization. Return how many builds (upper bound) have to happen next turn.
Definition safelanes.c:311
static const double LAMBDA
Definition safelanes.c:51
uint32_t FactionMask
A set of lane-building factions, represented as a bitfield.
Definition safelanes.c:89
static int * lane_faction
Definition safelanes.c:102
static const double JUMP_CONDUCTIVITY
Definition safelanes.c:52
void safelanes_recalculate(void)
Update the safe lane locations in response to the universe changing (e.g., diff applied).
Definition safelanes.c:252
static int * sys_to_first_edge
Definition safelanes.c:100
static void safelanes_initStacks_faction(void)
Sets up the local stacks with entry per faction.
Definition safelanes.c:446
static cholmod_common C
Definition safelanes.c:95
static void safelanes_initOptimizer(void)
Initializes resources used by lane optimization.
Definition safelanes.c:284
static void safelanes_initStacks_vertex(void)
Sets up the local stacks with entry per vertex (or per jump).
Definition safelanes.c:353
static cholmod_sparse * QtQ
Definition safelanes.c:111
static void safelanes_initFTilde(void)
Sets up the fluxes matrix f~.
Definition safelanes.c:628
SafeLane * safelanes_get(int faction, int standing, const StarSystem *system)
Gets a set of safelanes for a faction and system.
Definition safelanes.c:190
static cholmod_dense * ncholmod_ddmult(cholmod_dense *A, int transA, cholmod_dense *B)
Dense times dense matrix. Return A*B, or A'*B if transA is true.
Definition safelanes.c:956
static int * sys_to_first_vertex
Definition safelanes.c:98
static FactionMask MASK_ONE_FACTION(int id)
A mask giving this faction (NOT faction_stack index) exclusive rights to build, if it's a lane-buildi...
Definition safelanes.c:905
static double * cmp_key_ref
Definition safelanes.c:115
static int safelanes_triangleTooFlat(const vec2 *m, const vec2 *n, const vec2 *p, double lmn)
Return true if this triangle is so flat that lanes from point m to point n aren't allowed.
Definition safelanes.c:847
static double * tmp_edge_conduct
Definition safelanes.c:107
int Edge[2]
An edge is a pair of vertex indices.
Definition safelanes.c:79
static int * tmp_anchor_vertices
Definition safelanes.c:108
static Vertex * vertex_stack
Definition safelanes.c:96
static void safelanes_destroyOptimizer(void)
Frees resources used by lane optimization.
Definition safelanes.c:296
static Edge * tmp_jump_edges
Definition safelanes.c:106
static int * tmp_spob_indices
Definition safelanes.c:105
static cholmod_dense * utilde
Definition safelanes.c:113
VertexType
Object type: like SafeLaneLocType, but with Naev stack indexing in mind.
Definition safelanes.c:69
static void safelanes_initQtQ(void)
Sets up the (Q*)Q matrix.
Definition safelanes.c:602
static int FACTION_ID_TO_INDEX(int id)
Return the faction_stack index corresponding to a faction ID, or -1.
Definition safelanes.c:890
StarSystem * system_getIndex(int id)
Get the system by its index.
Definition space.c:989
StarSystem * system_getAll(void)
Gets an array (array.h) of all star systems.
Definition space.c:879
double system_getPresence(const StarSystem *sys, int faction)
Get the presence of a faction in a system.
Definition space.c:4181
Description of a lane-building faction.
Definition faction.c:54
double lane_base_cost
Definition faction.c:89
double lane_length_per_presence
Definition faction.c:88
int id
Definition safelanes.c:83
int devmode
Definition conf.h:157
Describes a safe lane, patrolled by a faction, within a system.
Definition safelanes.h:24
int point_id[2]
Definition safelanes.h:27
int faction
Definition safelanes.h:25
SafeLaneLocType point_type[2]
Definition safelanes.h:26
double bonus
Definition space.h:69
double base
Definition space.h:68
int faction
Definition space.h:67
Represents a Space Object (SPOB), including and not limited to planets, stations, wormholes,...
Definition space.h:89
SpobPresence presence
Definition space.h:104
Disjoint set forest on {0, .., n-1}.
Definition union_find.h:7
Reference to a spob or jump point.
Definition object.c:52
int index
Definition safelanes.c:75
int system
Definition safelanes.c:73
VertexType type
Definition safelanes.c:74
Represents a 2d vector.
Definition vec2.h:32
double y
Definition vec2.h:34
double x
Definition vec2.h:33
int * unionfind_findall(UnionFind *uf)
Returns a designated representative of each subset in an array (array.h).
Definition union_find.c:56
int unionfind_find(UnionFind *uf, int x)
Finds the designated representative of the subset containing x.
Definition union_find.c:48
void unionfind_init(UnionFind *uf, int n)
Creates a UnionFind structure on {0, ..., n}.
Definition union_find.c:14
void unionfind_union(UnionFind *uf, int x, int y)
Declares x and y to be in the same subset.
Definition union_find.c:34
void unionfind_free(UnionFind *uf)
Frees resources associated with uf.
Definition union_find.c:25