PIPS
unaires.c File Reference
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include "linear_assert.h"
#include "boolean.h"
#include "arithmetique.h"
#include "vecteur.h"
+ Include dependency graph for unaires.c:

Go to the source code of this file.

Macros

#define MALLOC(s, t, f)   malloc(s)
 package vecteur - operations unaires More...
 
#define FREE(p, t, f)   free(p)
 

Functions

void vect_normalize (Pvecteur v)
 void vect_normalize(Pvecteur v): division de tous les coefficients de v par leur pgcd; "normalisation" de vecteur directeur a coefficient entier More...
 
void vect_add_elem (Pvecteur *pvect, Variable var, Value val)
 void vect_add_elem(Pvecteur * pvect, Variable var, Value val): addition d'un vecteur colineaire au vecteur de base var au vecteur vect More...
 
void vect_erase_var (Pvecteur *ppv, Variable v)
 void vect_erase_var(Pvecteur * ppv, Variable v): projection du vecteur *ppv selon la direction v (i.e. More...
 
void vect_chg_coeff (Pvecteur *ppv, Variable var, Value val)
 void vect_chg_coeff(Pvecteur *ppv, Variable var, Value val): mise de la coordonnee var du vecteur *ppv a la valeur val More...
 
void vect_chg_var (Pvecteur *ppv, Variable v_old, Variable v_new)
 void vect_chg_var(Pvecteur *ppv, Variable v_old, Variable v_new) replace the variable v_old by v_new More...
 
Variable vect_one_coeff_if_any (Pvecteur v)
 
Pvecteur vect_del_var (Pvecteur v_in, Variable var)
 Pvecteur vect_del_var(Pvecteur v_in, Variable var): allocation d'un nouveau vecteur egal a la projection de v_in selon la direction var (i.e. More...
 
Value vect_coeff (Variable var, Pvecteur vect)
 Variable vect_coeff(Variable var, Pvecteur vect): coefficient de coordonnee var du vecteur vect —> Soit evar le vecteur de base de nom var: More...
 
Value vect_coeff_sum (Pvecteur vect)
 Value vect_coeff_sum(Pvecteur vect): coefficient sum de tout les val de ce vecteur (devrait etre dans reduction? FC) More...
 
Pvecteur vect_sign (Pvecteur v)
 Pvecteur vect_sign(Pvecteur v): application de l'operation signe au vecteur v. More...
 
void vect_sort_in_place (Pvecteur *pv, int *compare)
 void vect_sort_in_place(pv, compare) Pvecteur *pv; int (*compare)(Pvecteur *, Pvecteur *); More...
 
Pvecteur vect_sort (Pvecteur v, int *compare)
 Pvecteur vect_sort(v, compare) Pvecteur v; int (*compare)();. More...
 
int vect_compare (Pvecteur *pv1, Pvecteur *pv2)
 for qsort, returns: More...
 
void Pvecteur_separate_on_sign (Pvecteur v, Pvecteur *pvpos, Pvecteur *pvneg)
 void Pvecteur_separate_on_sign(v, pvpos, pvneg) Pvecteur v, *pvpos, *pvneg; More...
 
bool vect_common_variables_p (Pvecteur v1, Pvecteur v2)
 bool vect_common_variables_p(Pvecteur v1, v2) BA 19/05/94 input : two vectors. More...
 
bool vect_contains_variable_p (Pvecteur v, Variable var)
 bool vect_contains_variable_p(Pvecteur v, Variable var) BA 19/05/94 input : a vector and a variable output : true if var appears as a component of v, false otherwise. More...
 
int vect_lexicographic_compare (Pvecteur v1, Pvecteur v2, int(*compare)(Pvecteur *, Pvecteur *))
 qsort() is not safe if the comparison function is not antisymmetric. More...
 
int vect_lexicographic_compare2 (Pvecteur v1, Pvecteur v2, int(*compare)(Pvecteur *, Pvecteur *))
 Version for inequalities. More...
 
int vect_lexicographic_unsafe_compare (Pvecteur v1, Pvecteur v2, int(*compare)(Pvecteur *, Pvecteur *))
 
int vect_lexicographic_unsafe_compare2 (Pvecteur v1, Pvecteur v2, int(*compare)(Pvecteur *, Pvecteur *))
 
static int vect_lexicographic_coefficient_comparison (Pvecteur v1, Pvecteur v2)
 The two sparse vectors are assumed to have exactly the same structure, the same non-zero components in the same order. More...
 
int vect_lexicographic_unsafe_compare_generic (Pvecteur v1, Pvecteur v2, int(*compare)(Pvecteur *, Pvecteur *), bool is_equation __attribute__((unused)))
 This function is a trade-off between a real lexicographic sort and a prettyprinting function also used for code generation. More...
 

Macro Definition Documentation

◆ FREE

#define FREE (   p,
  t,
  f 
)    free(p)

Definition at line 42 of file unaires.c.

◆ MALLOC

#define MALLOC (   s,
  t,
  f 
)    malloc(s)

package vecteur - operations unaires

INTLIBRARY

Definition at line 41 of file unaires.c.

Function Documentation

◆ Pvecteur_separate_on_sign()

void Pvecteur_separate_on_sign ( Pvecteur  v,
Pvecteur pvpos,
Pvecteur pvneg 
)

void Pvecteur_separate_on_sign(v, pvpos, pvneg) Pvecteur v, *pvpos, *pvneg;

IN: v

OUT: pvpos, pvneg

this function builds 2 vectors composed of the positive and negative parts of the initial vector v which is not modified.

(c) FC 16/05/94

Parameters
pvposvpos
pvnegvneg

Definition at line 369 of file unaires.c.

371 {
372  Pvecteur vc;
373  Value val;
374  Variable var;
375 
376  *pvneg = VECTEUR_NUL,
377  *pvpos = VECTEUR_NUL;
378 
379  for(vc=v; vc; vc=vc->succ)
380  {
381  var = var_of(vc),
382  val = val_of(vc);
383  if (value_neg_p(val))
384  vect_add_elem(pvneg, var, value_uminus(val));
385  else
386  vect_add_elem(pvpos, var, val);
387  }
388 }
#define value_uminus(val)
unary operators on values
int Value
#define value_neg_p(val)
le type des coefficients dans les vecteurs: Value est defini dans le package arithmetique
Definition: vecteur-local.h:89
struct Svecteur * succ
Definition: vecteur-local.h:92
#define val_of(varval)
#define VECTEUR_NUL
DEFINITION DU VECTEUR NUL.
void * Variable
arithmetique is a requirement for vecteur, but I do not want to inforce it in all pips files....
Definition: vecteur-local.h:60
#define var_of(varval)
void vect_add_elem(Pvecteur *pvect, Variable var, Value val)
void vect_add_elem(Pvecteur * pvect, Variable var, Value val): addition d'un vecteur colineaire au ve...
Definition: unaires.c:72

References Svecteur::succ, val_of, value_neg_p, value_uminus, var_of, vect_add_elem(), and VECTEUR_NUL.

Referenced by Pcontrainte_to_expression_list().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ vect_add_elem()

void vect_add_elem ( Pvecteur pvect,
Variable  var,
Value  val 
)

void vect_add_elem(Pvecteur * pvect, Variable var, Value val): addition d'un vecteur colineaire au vecteur de base var au vecteur vect

--—> --—> —> *pvect := *pvect + val evar

Parameters
pvectvect
varar
valal

Definition at line 72 of file unaires.c.

73 {
74  if (val!=0)
75  {
76  Pvecteur vect;
77  for (vect=(*pvect); vect!=NULL; vect=vect->succ)
78  {
79  if (var_of(vect)==var)
80  {
81  value_addto(val_of(vect), val);
82  if (value_zero_p(val_of(vect)))
83  vect_erase_var(pvect, var_of(vect));
84  return;
85  }
86  }
87  // else: le coefficient valait 0 et n'etait pas represente
88  *pvect = vect_chain(*pvect, var, val);
89  }
90  // sinon, le vecteur est inchange et on ne fait rien
91 }
#define value_zero_p(val)
#define value_addto(ref, val)
Pvecteur vect_chain(Pvecteur v_in, Variable var, Value coeff)
package vecteur routines internes au package
Definition: private.c:69
void vect_erase_var(Pvecteur *ppv, Variable v)
void vect_erase_var(Pvecteur * ppv, Variable v): projection du vecteur *ppv selon la direction v (i....
Definition: unaires.c:106

References Svecteur::succ, val_of, value_addto, value_zero_p, var_of, vect_chain(), and vect_erase_var().

Referenced by add_affine_bound_conditions(), add_declaration_list_information(), add_elem_to_list_of_Pvecteur(), add_equivalence_equality(), add_fin_mat(), add_loop_index_exit_value(), add_loop_skip_condition(), add_var_sup(), adg_list_to_vect(), affine_to_transformer(), array_indices_communication(), base_reversal(), bitwise_xor_to_transformer(), bounds_equal_p(), build_image_base(), build_sc_machine(), cell_reference_sc_exact_projection_along_variable(), check_range_wrt_precondition(), complex_bound_computation(), compute_x_and_y_bounds(), constraints_to_loop_bound(), contrainte_reverse(), converti_psysmin_psysmax(), cout_nul(), dependence_cone_positive(), dependence_system_add_lci_and_di(), do_solve_hardware_constraints_on_nb_proc(), ecrit_coeff1(), eval_var(), expression_flt(), expression_less_than_in_context(), expression_multiply_sizeof_to_transformer(), find_motif(), find_pattern(), fonct_max_all(), fonct_max_d(), fonct_min_all(), fonct_min_d(), fonct_read(), formal_and_actual_parameters_association(), free_guards(), full_linearization(), gcd_and_constant_dependence_test(), generic_abs_to_transformer(), generic_equality_to_transformer(), generic_minmax_to_transformer(), hpfc_compute_lid(), iabs_to_transformer(), integer_divide_to_transformer(), integer_left_shift_to_transformer(), integer_minmax_to_transformer(), integer_multiply_to_transformer(), integer_power_to_transformer(), integer_right_shift_to_transformer(), invariant_vector_p(), list_to_base(), local_tile_constraints(), logical_binary_function_to_transformer(), logical_binary_operation_to_transformer(), logical_constant_to_transformer(), logical_unary_operation_to_transformer(), loop_bound_evaluation_to_transformer(), loop_bounds_to_tile_bounds(), loop_index_domaine_to_contrainte(), loop_regions_normalize(), lower_bound_generation(), make_base_of_nest(), make_constraint_expression(), make_context_of_loop(), make_datum_movement(), make_loop_indice_equation(), make_movements_loop_body_wp65(), make_scanning_over_one_tile(), make_tile_constraints(), my_system_remove_variables(), my_vect_substract(), new_constraint_for_coefficient_reduction_with_bounding_box(), partial_linearization(), pivoter_pas(), plint_degen(), polynome_roots(), polynome_sscanf(), polynome_to_vecteur(), Pvecteur_separate_on_sign(), reduce_loop_bound(), region_exact_projection_along_parameters(), region_exact_projection_along_variable(), relation_to_transformer(), remove_temporal_variables_from_system(), sc_add_di(), sc_add_dsi(), sc_bounded_normalization(), sc_find_equalities(), sc_minmax_of_variable_optim(), sc_multiply_constant_terms(), sc_of_constrs(), sc_proj_optim_on_di_ofl(), sc_projection_ofl_along_list_of_variables(), scanning_base_to_vect(), set_dimensions_of_local_variable_family(), simple_addition_to_transformer(), simplify_constraint_with_bounding_box(), simplify_float_constraint(), simplify_minmax_contrainte(), small_positive_slope_reduce_coefficients_with_bounding_box(), sys_int_redond(), test_bound_generation(), tile_change_of_basis(), tile_hyperplane_constraints(), tile_membership(), tile_membership_constraints(), transformer_add_3d_affine_constraint(), transformer_add_condition_information_updown(), transformer_add_equality(), transformer_add_equality_with_affine_term(), transformer_add_equality_with_integer_constant(), transformer_add_identity(), transformer_add_inequality(), transformer_add_inequality_with_affine_term(), transformer_add_inequality_with_integer_constraint(), transformer_add_integer_relation_information(), transformer_add_loop_index_initialization(), transformer_add_sign_information(), transformer_convex_hulls(), transformer_equality_fix_point(), transformer_logical_inequalities_add(), translate_global_value(), translate_to_module_frame(), update_lower_and_upper_bounds(), update_lower_or_upper_bound(), var_ecart_sup(), variables_in_declaration_list(), vars_read_and_written(), vect_add(), vect_cl_ofl_ctrl(), vect_gen_read(), vect_make(), vect_make_1D(), vect_make_dense(), vect_make_line(), vect_printout_order(), vect_reversal(), vect_substract(), vecteur_of_zvec(), vvs_to_sc(), xml_Chain_Graph(), xml_GlobalVariables(), xml_LocalVariables(), xml_Loops(), xml_Pattern_Paving(), xml_TaskParameters(), and xml_tiling().

+ Here is the call graph for this function:

◆ vect_chg_coeff()

void vect_chg_coeff ( Pvecteur ppv,
Variable  var,
Value  val 
)

void vect_chg_coeff(Pvecteur *ppv, Variable var, Value val): mise de la coordonnee var du vecteur *ppv a la valeur val

—> —> —> —> —> —> *ppv = *ppv - <*ppv . evar> evar + val evar

on n'a pas trouve de composante var

Parameters
ppvpv
varar
valal

Definition at line 143 of file unaires.c.

147 {
148  Pvecteur pvcour;
149 
150  if (val == 0) {
151  vect_erase_var(ppv, var);
152  }
153  else {
154  for (pvcour = (*ppv); pvcour != NULL; pvcour = pvcour->succ) {
155  if (pvcour->var == var) {
156  pvcour->val = val;
157  return;
158  }
159  }
160  /* on n'a pas trouve de composante var */
161  *ppv = vect_chain(*ppv,var,val);
162  }
163 }
Value val
Definition: vecteur-local.h:91
Variable var
Definition: vecteur-local.h:90

References Svecteur::succ, Svecteur::val, Svecteur::var, vect_chain(), and vect_erase_var().

Referenced by build_sc_with_several_uniform_ref(), complex_bound_computation(), contrainte_normalize(), cout_nul(), eval_var(), find_eg(), find_vbase(), free_guards(), get_exp_schedule(), gomory_trait_eq(), loop_nest_to_offset(), loop_nest_to_tile(), loop_nest_to_wp65_code(), lower_bound_generation(), mat_sys_conv(), matrice_index_sys(), matrices_to_constraints(), matrices_to_constraints_with_sym_cst(), matrices_to_contraintes_with_sym_cst(), matrices_to_loop_sc(), matrices_to_sc(), movement_computation(), my_contrainte_normalize(), my_matrices_to_constraints_with_sym_cst(), my_matrices_to_constraints_with_sym_cst_2(), oter_lvbase(), pivoter_pas(), plint_degen(), print_call_precondition(), pu_matrices_to_contraintes(), sc_consistent_p(), sc_image_computation(), sc_weak_consistent_p(), test_bound_generation(), transformer_add_variable_incrementation(), transformer_convex_hulls(), update_basis(), update_coefficient_signs_in_vector(), upper_bound_generation(), var_pivotd(), var_pivots(), var_posit(), vect_subst(), xml_Chain_Graph(), and xml_Compute_and_Need().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ vect_chg_var()

void vect_chg_var ( Pvecteur ppv,
Variable  v_old,
Variable  v_new 
)

void vect_chg_var(Pvecteur *ppv, Variable v_old, Variable v_new) replace the variable v_old by v_new

Parameters
ppvpv
v_old_old
v_new_new

Definition at line 168 of file unaires.c.

171 {
172  Pvecteur pvcour;
173 
174  for (pvcour = (*ppv); pvcour != NULL; pvcour = pvcour->succ) {
175  if (pvcour->var == v_old){
176  pvcour->var = v_new;
177  return;
178  }
179  }
180 }

References Svecteur::succ, and Svecteur::var.

Referenced by gcd_and_constant_dependence_test(), poly_chg_var(), polynome_chg_var(), polynome_roots(), and sc_chg_var().

+ Here is the caller graph for this function:

◆ vect_coeff()

Value vect_coeff ( Variable  var,
Pvecteur  vect 
)

Variable vect_coeff(Variable var, Pvecteur vect): coefficient de coordonnee var du vecteur vect —> Soit evar le vecteur de base de nom var:

    --->   --->

return <vect . evar>; (i.e. return vect[var])

Parameters
varar
vectect

Definition at line 228 of file unaires.c.

231 {
232  for ( ; vect != NULL ; vect = vect->succ)
233  if (var_of(vect) == var) {
234  assert(val_of(vect)!=VALUE_ZERO);
235  return(val_of(vect));
236  }
237  return VALUE_ZERO;
238 }
#define VALUE_ZERO
#define assert(ex)
Definition: newgen_assert.h:41

References assert, Svecteur::succ, val_of, VALUE_ZERO, and var_of.

Referenced by add_bounding_box_constraints(), add_declaration_list_information(), add_reference_information(), add_var_sup(), affine_expression_of_loop_index_p(), aligned_p(), alignment_p(), array_access_to_array_ranges(), array_indices_communication(), array_overflow(), atomize_one_message(), bound_distribution(), bound_redund_with_sc_p(), bounds_equal_p(), build_integer_sc_nredund(), build_list_of_min(), build_sc_machine(), build_third_comb(), calculate_delay(), call_rwt(), change_base_in_sc(), check_coefficient_reduction(), combiner_ofl_with_test(), compare_the_constraints(), complex_bound_computation(), compose_vvs(), compute_receive_content(), compute_receive_domain(), compute_x_and_y_bounds(), const_negative(), constraint_distribution(), constraint_integer_combination(), constraint_to_bound(), constraint_without_vars(), constraints_eliminate_constant_terms(), constraints_for_bounds(), constraints_to_loop_bound(), constraints_to_matrices(), constraints_with_sym_cst_to_matrices(), constrs_of_sc(), contrainte_dup_extract(), contrainte_eval(), contrainte_normalize(), contrainte_parallele(), contrainte_subst_ofl_ctrl(), contrainte_to_matrix_ligne(), contrainte_var_min_coeff(), contraintes_to_expression(), contraintes_with_sym_cst_to_matrices(), convert_bound_expression(), convex_in_effect_loop_range_fix(), cost_of_constant_operations(), cst_vector_p(), dj_variable_substitution_with_eqs_ofl_ctrl(), do_group_statement_constant(), ecrit_ligne(), eligible_for_coefficient_reduction_with_bounding_box_p(), elim_var_with_eg(), eq_diff_const(), eq_sum_const(), eq_var_nophi_min_coeff(), eq_var_phi(), eval(), eval_2D_vecteur(), eval_var(), expression_and_precondition_to_integer_interval(), expression_integer_constant_p(), expression_to_int(), find_eg(), find_first_integer_point_in_between(), find_motif(), find_pattern(), find_vbase(), free_guards(), gcd_and_constant_dependence_test(), generate_one_message(), get_bounds_expression(), get_const_diff(), get_const_off(), gomory_eq(), gomory_trait_eq(), hpfc_broadcast_buffers(), hpfc_compute_lid(), hpfc_integer_constant_expression_p(), HpfcExpressionToInt(), include_trans_on_LC_in_ref(), inegalite_comb_ofl_ctrl(), initialize_offsets(), intersection(), is_good_direction_p(), is_inferior_monome(), Lbound(), legal_point_p(), ligne_pivot(), lignes_entrant(), loop_bounds_to_tile_bounds(), loop_flt(), loop_index_in_several_indices(), loop_regions_normalize(), loop_sc_to_matrices(), lower_bound_generation(), make_bounds(), make_constraint_expression(), make_loop_indice_equation(), make_vvs_from_sc(), message_larger_p(), message_manageable_p(), monome_sprint(), my_constraints_with_sym_cst_to_matrices(), my_contrainte_normalize(), my_substitute_var_with_vec(), my_vect_var_subst(), new_constraint_for_coefficient_reduction_with_bounding_box(), new_ecrit_ligne(), one_message_guards_and_neighbour(), one_receive_message(), outliner_smart_references_computation(), partial_broadcast_coefficients(), phi_free_contraints_to_expressions(), pivoter(), pivoter_pas(), plc_elim_var_with_eg(), plc_make_dim(), plc_make_vvs_with_vector(), plint(), plint_pas(), polynome_contains_var(), polynome_degree(), polynome_factorize(), polynome_roots(), polynome_var_subst(), print_call_precondition(), print_cone_vecteur(), print_vect_in_vertice_val(), pu_contraintes_to_matrices(), Pvecteur_to_assign_statement(), pvecteur_to_polynome(), region_projection_along_index_safe_p(), region_range_nul_p(), sc_add_di(), sc_add_dsi(), sc_bounded_normalization(), sc_constrains_variable_p(), sc_elim_double_constraints(), sc_elim_triang_integer_redund_constraint_p(), sc_elim_var(), sc_eliminate_constant_terms(), sc_find_equalities(), sc_integer_projection_information(), sc_minmax_of_variable(), sc_minmax_of_variable_optim(), sc_normalize2(), sc_safe_elim_db_constraints(), sc_simplex_feasibility_ofl_ctrl_fixprec(), sc_to_iproblem(), sc_to_matrices(), sc_variable_rename(), separate_variables(), set_dimensions_of_local_variable_family(), sg_fprint_as_ddv(), shift_expression_of_loop_index_p(), signed_integer_constant_expression_value(), simple_indices_p(), simplify_constraint_with_bounding_box(), simplify_dimension(), simplify_minmax_contrainte(), simplify_predicate(), small_positive_slope_reduce_coefficients_with_bounding_box(), small_slope_and_first_quadrant_p(), sol_entiere(), sol_finale(), sol_positive(), sol_positive_simpl(), sort_tile_indices(), st_one_message(), stmt_bdt_directions(), substitute_var_with_vec(), supported_ref_p(), sys_mat_conv(), sys_matrice_index(), system_new_var_subst(), test_borne(), test_bound_generation(), tile_change_of_basis(), tile_hyperplane_constraints(), tile_membership_constraints(), transform_in_ineq(), transformer_affect_linear_p(), translate_to_module_frame(), Ubound(), update_coefficient_signs_in_vector(), update_indices_for_local_computation(), update_lower_and_upper_bounds(), update_lower_or_upper_bound(), upper_bound_generation(), var_in_lcontrainte_p(), var_pivotd(), var_pivots(), var_with_unity_coeff_p(), vars_in_vect_p(), vect_const_p(), vect_equal(), vect_equal_except(), vect_fprint_as_dense(), vect_in_p(), vect_oppos(), vect_opposite_except(), vect_prod_scal(), vect_product(), vect_proport(), vect_sprint_as_monome(), vect_subst(), vect_substitute_dimension(), vect_var_subst(), vvs_on_vecteur(), vvs_on_vvs(), which_array_dimension(), xml_Chain_Graph(), xml_Compute_and_Need(), xml_ConstOffset(), xml_GlobalVariables(), xml_LocalVariables(), xml_LoopOffset(), xml_Pattern_Paving(), xml_Region_Range(), xml_TaskParameter(), xml_tiling(), and zmat_set_row().

◆ vect_coeff_sum()

Value vect_coeff_sum ( Pvecteur  vect)

Value vect_coeff_sum(Pvecteur vect): coefficient sum de tout les val de ce vecteur (devrait etre dans reduction? FC)

return Value Lei Zhou Mar.25, 91

Parameters
vectect

Definition at line 246 of file unaires.c.

248 {
249  Value val = VALUE_ZERO;
250 
251  if ( vect->var == TCST )
252  return val;
253  for (; vect != NULL ; vect = vect->succ) {
254  value_addto(val,vecteur_val(vect));
255  assert(value_notzero_p(val_of(vect)));
256  }
257  return val;
258 }
#define value_notzero_p(val)
#define TCST
VARIABLE REPRESENTANT LE TERME CONSTANT.
#define vecteur_val(v)

References assert, Svecteur::succ, TCST, val_of, value_addto, value_notzero_p, VALUE_ZERO, Svecteur::var, and vecteur_val.

Referenced by is_inferior_monome().

+ Here is the caller graph for this function:

◆ vect_common_variables_p()

bool vect_common_variables_p ( Pvecteur  v1,
Pvecteur  v2 
)

bool vect_common_variables_p(Pvecteur v1, v2) BA 19/05/94 input : two vectors.

output : true if they have at least one common variable, false otherwise. modifies : nothing.

Parameters
v11
v22

Definition at line 397 of file unaires.c.

399 {
400  Pvecteur ev;
401 
402  for(ev = v1; !VECTEUR_NUL_P(ev); ev = ev->succ) {
404  return true;
405  }
406  return false;
407 }
#define vecteur_var(v)
#define VECTEUR_NUL_P(v)
bool vect_contains_variable_p(Pvecteur v, Variable var)
bool vect_contains_variable_p(Pvecteur v, Variable var) BA 19/05/94 input : a vector and a variable o...
Definition: unaires.c:415

References Svecteur::succ, vect_contains_variable_p(), VECTEUR_NUL_P, and vecteur_var.

Referenced by elementary_convex_union(), and sc_restricted_to_variables_transitive_closure().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ vect_compare()

int vect_compare ( Pvecteur pv1,
Pvecteur pv2 
)

for qsort, returns:

  • if v1 < v2 0 if v1 = v2
  • if v1 > v2
Parameters
pv1v1
pv2v2

Definition at line 352 of file unaires.c.

354 {
355  return(strcmp((char *)&var_of(*pv1), (char *)&var_of(*pv2)));
356 }

References var_of.

Referenced by find_vbase(), MaxBoundary(), MinBoundary(), polynome_sort(), var_pivotd(), and var_pivots().

+ Here is the caller graph for this function:

◆ vect_contains_variable_p()

bool vect_contains_variable_p ( Pvecteur  v,
Variable  var 
)

bool vect_contains_variable_p(Pvecteur v, Variable var) BA 19/05/94 input : a vector and a variable output : true if var appears as a component of v, false otherwise.

modifies : nothing

Parameters
varar

Definition at line 415 of file unaires.c.

418 {
419  bool in_base;
420 
421  for(; !VECTEUR_NUL_P(v) && !variable_equal(vecteur_var(v), var); v = v->succ)
422  ;
423  in_base = !VECTEUR_NUL_P(v);
424  return(in_base);
425 }
bool variable_equal(Variable v1, Variable v2)
package vecteur - routines sur les variables
Definition: variable.c:62

References variable_equal(), VECTEUR_NUL_P, and vecteur_var.

Referenced by DivExists(), translate_to_module_frame(), vect_common_variables_p(), and vect_same_variables_p().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ vect_del_var()

Pvecteur vect_del_var ( Pvecteur  v_in,
Variable  var 
)

Pvecteur vect_del_var(Pvecteur v_in, Variable var): allocation d'un nouveau vecteur egal a la projection de v_in selon la direction var (i.e.

le coefficient de la coordonnee var est mis a 0)

Soit evar le vecteur de base correspondant a var:

      ---->

allocate v_out;

-—> —> -—> —> —> v_out := v_in - <v_out . evar> evar

  ---->

return v_out;

Parameters
v_in_in
varar

Definition at line 206 of file unaires.c.

209 {
210  if(v_in!=NULL){
211  Pvecteur v_out = vect_dup(v_in);
212  vect_erase_var(&v_out,var);
213  return(v_out);
214  }
215  else
216  return(NULL);
217 }
Pvecteur vect_dup(Pvecteur v_in)
Pvecteur vect_dup(Pvecteur v_in): duplication du vecteur v_in; allocation de et copie dans v_out;.
Definition: alloc.c:51

References vect_dup(), and vect_erase_var().

Referenced by affine_expression_of_loop_index_p(), array_access_to_array_ranges(), build_list_of_min(), constraint_to_bound(), constraints_to_loop_bound(), contraintes_to_expression(), do_group_statement_constant(), get_bounds_expression(), hpfc_integer_constant_expression_p(), make_vvs_from_sc(), monome_del_var(), outliner_smart_references_computation(), polynome_roots(), shift_expression_of_loop_index_p(), the_index_of_vect(), and translate_to_module_frame().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ vect_erase_var()

void vect_erase_var ( Pvecteur ppv,
Variable  v 
)

void vect_erase_var(Pvecteur * ppv, Variable v): projection du vecteur *ppv selon la direction v (i.e.

mise a zero de la coordonnee v du vecteur pointe par ppv)

Soit ev le vecteur de base correspondant a v:

—> —> —> -> -> *ppv := *ppv - <*ppv . ev> ev

Note: cette routine ne fait pas l'hypothese que chaque coordonnee n'apparait qu'une fois; on pourrait l'accelerer en forcant pvcour a NULL des que la coordonnee est trouvee.

A-t-on trouve la composante v?

Si oui, est-il possible de la dechainer?

elle n'est pas en tete de liste

Elle est en tete de liste; il faut modifier ppv

Non, on passe a la composante suivante...

Parameters
ppvpv

Definition at line 106 of file unaires.c.

109 {
110  Pvecteur pvprec, pvcour;
111 
112  for (pvprec = NULL, pvcour = (*ppv); pvcour != NULL;) {
113  /* A-t-on trouve la composante v? */
114  if (pvcour->var == v) {
115  /* Si oui, est-il possible de la dechainer? */
116  if (pvprec != NULL) {
117  /* elle n'est pas en tete de liste */
118  Pvecteur pvprim = pvcour;
119  pvcour = pvprec->succ = pvcour->succ;
120  FREE((char *)pvprim, VECTEUR, "vect_erase_var");
121  }
122  else {
123  /* Elle est en tete de liste; il faut modifier ppv */
124  *ppv = pvcour->succ;
125  FREE((char *)pvcour,VECTEUR,"vect_erase_var");
126  pvcour = *ppv;
127  }
128  }
129  else {
130  /* Non, on passe a la composante suivante... */
131  pvprec = pvcour;
132  pvcour = pvcour->succ;
133  }
134  }
135 }
#define VECTEUR
package sur les vecteurs creux et les bases
Definition: vecteur-local.h:51
#define FREE(p, t, f)
Definition: unaires.c:42

References FREE, Svecteur::succ, Svecteur::var, and VECTEUR.

Referenced by algorithm_row_echelon_generic(), base_remove_variable(), build_third_comb(), change_base_in_sc(), compose_vvs(), constraints_to_loop_bound(), find_motif(), find_pattern(), free_guards(), include_trans_on_LC_in_ref(), is_inferior_monome(), make_constraint_expression(), movement_computation(), my_substitute_var_with_vec(), my_vect_var_subst(), new_system_with_only_live_variable(), pip_solve(), pip_solve_min_with_big(), Pvecteur_to_assign_statement(), sc_find_equalities(), sc_fm_project_variables(), sc_force_variable_to_zero(), simplify_minmax_contrainte(), system_new_var_subst(), the_index_of_vect(), transform_in_ineq(), update_basis(), vect_add_elem(), vect_add_first(), vect_change_base(), vect_chg_coeff(), vect_del_var(), vect_printout_order(), vect_var_subst(), vvs_on_vecteur(), vvs_on_vvs(), xml_Boxes(), and xml_Region_Range().

+ Here is the caller graph for this function:

◆ vect_lexicographic_coefficient_comparison()

static int vect_lexicographic_coefficient_comparison ( Pvecteur  v1,
Pvecteur  v2 
)
static

The two sparse vectors are assumed to have exactly the same structure, the same non-zero components in the same order.

Definition at line 486 of file unaires.c.

487 {
488  Pvecteur pv1, pv2;
489  int cmp = 0;
490 
491  for(pv1 = v1, pv2 = v2;
492  !VECTEUR_UNDEFINED_P(pv1) && !VECTEUR_UNDEFINED_P(pv2) && cmp == 0;
493  pv1 = pv1->succ, pv2= pv2->succ) {
494  Value c1 = vecteur_val(pv1);
495  Value c2 = vecteur_val(pv2);
496  if(value_gt(c1, c2))
497  cmp = 1;
498  else if(value_gt(c2, c1))
499  cmp = -1;
500  }
501  return cmp;
502 }
#define value_gt(v1, v2)
#define VECTEUR_UNDEFINED_P(v)

References Svecteur::succ, value_gt, VECTEUR_UNDEFINED_P, and vecteur_val.

Referenced by vect_lexicographic_unsafe_compare_generic().

+ Here is the caller graph for this function:

◆ vect_lexicographic_compare()

int vect_lexicographic_compare ( Pvecteur  v1,
Pvecteur  v2,
int(*)(Pvecteur *, Pvecteur *)  compare 
)

qsort() is not safe if the comparison function is not antisymmetric.

It wanders out of the array to be sorted. It's a pain to debug. Let's play safe. Version for equations

Parameters
v11
v22

Definition at line 433 of file unaires.c.

435 {
436  int cmp12 = 0;
437  int cmp21 = 0;
438 
439  cmp12 = vect_lexicographic_unsafe_compare(v1, v2, compare);
440  cmp21 = vect_lexicographic_unsafe_compare(v2, v1, compare);
441 
442  assert(cmp12 == -cmp21);
443 
444  return cmp12;
445 }
int vect_lexicographic_unsafe_compare(Pvecteur v1, Pvecteur v2, int(*compare)(Pvecteur *, Pvecteur *))
Definition: unaires.c:464

References assert, and vect_lexicographic_unsafe_compare().

Referenced by equation_lexicographic_compare(), and sc_lexicographic_sort().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ vect_lexicographic_compare2()

int vect_lexicographic_compare2 ( Pvecteur  v1,
Pvecteur  v2,
int(*)(Pvecteur *, Pvecteur *)  compare 
)

Version for inequalities.

Parameters
v11
v22

Definition at line 449 of file unaires.c.

451 {
452  int cmp12 = 0;
453  int cmp21 = 0;
454 
455  cmp12 = vect_lexicographic_unsafe_compare2(v1, v2, compare);
456  cmp21 = vect_lexicographic_unsafe_compare2(v2, v1, compare);
457 
458  assert(cmp12 == -cmp21);
459 
460  return cmp12;
461 }
int vect_lexicographic_unsafe_compare2(Pvecteur v1, Pvecteur v2, int(*compare)(Pvecteur *, Pvecteur *))
Definition: unaires.c:474

References assert, and vect_lexicographic_unsafe_compare2().

Referenced by inequality_lexicographic_compare().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ vect_lexicographic_unsafe_compare()

int vect_lexicographic_unsafe_compare ( Pvecteur  v1,
Pvecteur  v2,
int(*)(Pvecteur *, Pvecteur *)  compare 
)
Parameters
v11
v22

Definition at line 464 of file unaires.c.

466 {
467  int cmp = 0;
468 
469  cmp = vect_lexicographic_unsafe_compare_generic(v1, v2, compare, true);
470 
471  return cmp;
472 }
int vect_lexicographic_unsafe_compare_generic(Pvecteur v1, Pvecteur v2, int(*compare)(Pvecteur *, Pvecteur *), bool is_equation __attribute__((unused)))
This function is a trade-off between a real lexicographic sort and a prettyprinting function also use...
Definition: unaires.c:542

References vect_lexicographic_unsafe_compare_generic().

Referenced by vect_lexicographic_compare().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ vect_lexicographic_unsafe_compare2()

int vect_lexicographic_unsafe_compare2 ( Pvecteur  v1,
Pvecteur  v2,
int(*)(Pvecteur *, Pvecteur *)  compare 
)
Parameters
v11
v22

Definition at line 474 of file unaires.c.

476 {
477  int cmp = 0;
478 
479  cmp = vect_lexicographic_unsafe_compare_generic(v1, v2, compare, false);
480 
481  return cmp;
482 }

References vect_lexicographic_unsafe_compare_generic().

Referenced by vect_lexicographic_compare2().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ vect_lexicographic_unsafe_compare_generic()

int vect_lexicographic_unsafe_compare_generic ( Pvecteur  v1,
Pvecteur  v2,
int(*)(Pvecteur *, Pvecteur *)  compare,
bool is_equation   __attribute__(unused) 
)

This function is a trade-off between a real lexicographic sort and a prettyprinting function also used for code generation.

Its goal is to sort constraints used as loop bounds or obtained as preconditions. Each constraint or vector is assumed internally sorted using the compare function, e.g. 2i + 3j + 4 is correct, while 3j + 4 + 2i is not if the alphabetical order is used and if constants appear as last vector elements.

When multiple constraints appear in a constraint system, we usually want simpler constraints and vectors first but some order between variables is still used. E.g. {1, i, i+1, i+j+1, j}. Here, we do not take into account the vector length as a primary critarion, but the alphabetical order. The above system is first reduced to a set of "words" {"", "i", "ij", "j"} and is lexicographically sorted [This is simplified: each letter in this example is in fact a work in the general case]. In case, two "words" are identical, e.g. "ij" and "ij", the length of the two underlying vectors, e.g. i+j, i+j+1, i+2*j, are compared. If they are equal, the lexicographic order of the coefficients is used to disambiguate the comparison, e.g. i+j < i+2j.

A lot of problems arise because 0 is not represented. It's hard to compare I==0 and I==1 because they do not have the same numbers of sparse components.

Furthermore, vectors representing equalities and inequalities must not be handled in the same way because only equalities can be multiplied by -1.

Not satisfying for Transformations/Tiling.sub/tiling05: 0 constant terms should be handled in a special way so as to have "i" be less than "i-1". Since -1 is less than 0, the longest constraint comes first.

Not satisfying for Transformations/Tiling.sub/tiling04: scopes had an impact on the comparisons that is not natural for users.

It is assumed that vectors v1 and v2 are already lexicographically sorted, according to the same lexicographic order.

The constant term should always be the last one. But the lexicographic sort now move them ahead! So we have to skip them at each position but we assume that they are both either leading or trailing elements.

Lexicographic comparison on variable names

Use vector lengths as discriminator in case an initial constant term makes a difference.

Use lexicographic order on coefficients as discriminator

This point should not be reachable unless the two vectors are identical and end with two constant terms

We need a total order to avoid non-determinism, i.e. dependence on pointer values

Definition at line 542 of file unaires.c.

545 {
546  /* It is assumed that vectors v1 and v2 are already
547  lexicographically sorted, according to the same lexicographic
548  order.
549 
550  The constant term should always be the last one. But the
551  lexicographic sort now move them ahead! So we have to skip them
552  at each position but we assume that they are both either
553  leading or trailing elements. */
554  int cmp = 0;
555  Pvecteur pv1 = term_cst(v1)? v1->succ : v1;
556  Pvecteur pv2 = term_cst(v2)? v2->succ : v2;
557 
558  /* Lexicographic comparison on variable names */
559  for(;
560  !VECTEUR_UNDEFINED_P(pv1) && !VECTEUR_UNDEFINED_P(pv2) && cmp == 0
561  && !term_cst(pv1) && !term_cst(pv2);
562  pv1 = pv1->succ, pv2= pv2->succ) {
563 
564  cmp = compare(&pv1, &pv2);
565  }
566 
567  if(cmp==0) {
568  if(VECTEUR_UNDEFINED_P(pv1)) {
569  if(VECTEUR_UNDEFINED_P(pv2)) {
570  /* Use vector lengths as discriminator in case an initial
571  constant term makes a difference. */
572  int n1 = vect_size(v1);
573  int n2 = vect_size(v2);
574  if(n1>n2)
575  cmp = 1;
576  else if(n1<n2)
577  cmp = -1;
578  else {
579  /* Use lexicographic order on coefficients as discriminator */
581  }
582  }
583  else
584  cmp = -1; // v2 is longer
585  }
586  else {
587  if(VECTEUR_UNDEFINED_P(pv2))
588  cmp = 1; // v1 is longer
589  else {
590  if(term_cst(pv1)) {
591  if(term_cst(pv2)) {
592  // Use only constant terms as differentiator
593  // cmp = value_comparison(vecteur_val(pv1), vecteur_val(pv2));
594  // To get lower bounds before upper bounds when only one variables is constrained
596  }
597  else
598  cmp = -1;
599  }
600  else {
601  if(term_cst(pv2))
602  cmp = 1;
603  else
604  /* This point should not be reachable unless the two vectors
605  are identical and end with two constant terms */
606  cmp = 0;
607  }
608  }
609  }
610  }
611 
612  /* We need a total order to avoid non-determinism, i.e. dependence
613  on pointer values */
614  assert(cmp!=0 || vect_equal(v1, v2));
615 
616  return cmp;
617 }
bool vect_equal(Pvecteur v1, Pvecteur v2)
bool vect_equal(Pvecteur v1, Pvecteur v2): test a egalite de deux vecteurs
Definition: reductions.c:278
int vect_size(Pvecteur v)
package vecteur - reductions
Definition: reductions.c:47
#define term_cst(varval)
static int vect_lexicographic_coefficient_comparison(Pvecteur v1, Pvecteur v2)
The two sparse vectors are assumed to have exactly the same structure, the same non-zero components i...
Definition: unaires.c:486

References assert, Svecteur::succ, term_cst, vect_equal(), vect_lexicographic_coefficient_comparison(), vect_size(), and VECTEUR_UNDEFINED_P.

Referenced by vect_lexicographic_unsafe_compare(), and vect_lexicographic_unsafe_compare2().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ vect_normalize()

void vect_normalize ( Pvecteur  v)

void vect_normalize(Pvecteur v): division de tous les coefficients de v par leur pgcd; "normalisation" de vecteur directeur a coefficient entier

unaires.c

-> -> -> -> si v == 0 alors v := 0; sinon pgcd = PGCD v[i]; i -> -> v := v / pgcd

Le pgcd est toujours positif.

Ancien nom: vect_norm()

Definition at line 59 of file unaires.c.

60 {
61  Value gcd = vect_pgcd_all(v);
62  if (value_notzero_p(gcd) && value_notone_p(gcd))
63  (void) vect_div(v, gcd);
64 }
#define value_notone_p(val)
Value vect_pgcd_all(Pvecteur v)
Value vect_pgcd(Pvecteur v): calcul du pgcd de tous les coefficients non nul d'un vecteur v.
Definition: reductions.c:108
Pvecteur vect_div(Pvecteur v, Value x)
Pvecteur vect_div(Pvecteur v, Value x): division du vecteur v par le scalaire x, si x est different d...
Definition: scalaires.c:52

References value_notone_p, value_notzero_p, vect_div(), and vect_pgcd_all().

Referenced by bounds_equal_p(), build_third_comb(), compute_x_and_y_bounds(), include_trans_on_LC_in_ref(), make_rational_exp(), my_sc_normalize(), norm_eq(), ray_dte_normalize(), sc_add_normalize_eq(), sc_add_normalize_ineq(), sc_elim_double_constraints(), sc_elim_redund(), sc_elim_redund_with_first_ofl_ctrl(), sc_gcd_normalize(), sc_normalize(), simplify_relational_expression(), system_new_var_subst(), and vect_make_line().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ vect_one_coeff_if_any()

Variable vect_one_coeff_if_any ( Pvecteur  v)

Definition at line 182 of file unaires.c.

183 {
184  for (; v; v=v->succ)
185  if (v->var && (value_one_p(v->val) || value_mone_p(v->val)))
186  return v->var;
187  return NULL;
188 }
#define value_mone_p(val)
#define value_one_p(val)

References Svecteur::succ, Svecteur::val, value_mone_p, value_one_p, and Svecteur::var.

◆ vect_sign()

Pvecteur vect_sign ( Pvecteur  v)

Pvecteur vect_sign(Pvecteur v): application de l'operation signe au vecteur v.

-> -> v := signe(v ); -> return v ;

Definition at line 269 of file unaires.c.

271 {
272  Pvecteur coord;
273 
274  for(coord = v; coord!=NULL; coord=coord->succ)
275  val_of(coord) = int_to_value(value_sign(val_of(coord)));
276 
277  return v;
278 }
#define value_sign(v)
trian operators on values
#define int_to_value(i)
end LINEAR_VALUE_IS_INT

References int_to_value, Svecteur::succ, val_of, and value_sign.

◆ vect_sort()

Pvecteur vect_sort ( Pvecteur  v,
int compare 
)

Pvecteur vect_sort(v, compare) Pvecteur v; int (*compare)();.

—> --> OUT = sorted IN

Definition at line 335 of file unaires.c.

338 {
339  Pvecteur
340  new = vect_dup(v);
341 
342  vect_sort_in_place(&new, compare);
343  return(new);
344 }
void vect_sort_in_place(Pvecteur *pv, int *compare)
void vect_sort_in_place(pv, compare) Pvecteur *pv; int (*compare)(Pvecteur *, Pvecteur *);
Definition: unaires.c:290

References vect_dup(), and vect_sort_in_place().

Referenced by constraints_to_loop_bound(), find_vbase(), is_inferior_monome(), make_vecteur_expression(), polynome_sort(), polynome_used_var(), var_pivotd(), and var_pivots().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ vect_sort_in_place()

void vect_sort_in_place ( Pvecteur pv,
int compare 
)

void vect_sort_in_place(pv, compare) Pvecteur *pv; int (*compare)(Pvecteur *, Pvecteur *);

Sorts the vector in place. It is an interface to qsort (stdlib). see man qsort about the compare function, which tells < == or >.

FC 29/12/94

the temporary table is created and initialized

sort!

FI: I do not know how to cast compare() properly

qsort(table, n, sizeof(Pvecteur), int (* compare)());

the vector is regenerated in order

clean and return

Definition at line 290 of file unaires.c.

293 {
294  int
295  n = vect_size(*pv);
296  Pvecteur
297  v,
298  *table,
299  *point;
300 
301  if ( n==0 || n==1 ) return;
302 
303  /* the temporary table is created and initialized
304  */
305  table = (Pvecteur*) malloc(sizeof(Pvecteur)*n);
306 
307  for (v=*pv, point=table; v!=(Pvecteur)NULL; v=v->succ, point++)
308  *point=v;
309 
310  /* sort!
311  */
312  /* FI: I do not know how to cast compare() properly */
313  /* qsort(table, n, sizeof(Pvecteur), int (* compare)()); */
314  qsort(table, n, sizeof(Pvecteur),(int (*)()) compare);
315 
316  /* the vector is regenerated in order
317  */
318  for (point=table; n>1; point++, n--)
319  (*point)->succ=*(point+1);
320 
321  (*point)->succ=(Pvecteur) NULL;
322 
323  /* clean and return
324  */
325  *pv=*table; free(table);
326 }
void * malloc(YYSIZE_T)
void free(void *)
struct Svecteur * Pvecteur

References free(), malloc(), Svecteur::succ, and vect_size().

Referenced by contrainte_vect_sort(), sc_find_equalities(), sc_lexicographic_sort(), transformer_normalize(), vect_printout_order(), and vect_sort().

+ Here is the call graph for this function:
+ Here is the caller graph for this function: