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

Go to the source code of this file.

Macros

#define MALLOC(s, t, f)   malloc(s)
 INTLIBRARY. More...
 
#define FREE(p, t, f)   free(p)
 

Functions

Pvecteur vect_dup (Pvecteur v_in)
 Pvecteur vect_dup(Pvecteur v_in): duplication du vecteur v_in; allocation de et copie dans v_out;. More...
 
void vect_rm (Pvecteur v)
 void vect_rm(Pvecteur v): desallocation des couples de v; More...
 
Pvecteur vect_new (Variable var, Value coeff)
 Pvecteur vect_new(Variable var,Value coeff): allocation d'un vecteur colineaire au vecteur de base var et de coefficient coeff (i.e. More...
 
void dbg_vect_rm (Pvecteur v, char __attribute__((unused)) *f)
 void dbg_vect_rm(Pvecteur v, char * f): desallocation d'un vecteur avec marquage de la fonction provoquant la desallocation More...
 
Pvecteur vect_make (Pvecteur v, Variable var, Value val,...)
 Pvecteur vect_make(v, [var, val,]* 0, val) Pvecteur v; // may be NULL, use assigne anyway Variable var; Value val;. More...
 
Pvecteur vect_make_dense (Pbase b, Value val,...)
 Allocate a new vector v whose coefficient are given by the list of values ad whose dimension is given by b. More...
 
Pvecteur vect_make_1D (Value a, Variable x, Value b)
 Generate a sparse vector a x + b TCST. More...
 
Pbase vect_copy (Pvecteur b)
 direct duplication. More...
 
Pbase base_dup (Pbase b)
 Pbase base_dup(Pbase b) Note: this function changes the value of the pointer. More...
 
Pbase base_copy (Pbase b)
 Direct duplication. More...
 

Macro Definition Documentation

◆ FREE

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

Definition at line 43 of file alloc.c.

◆ MALLOC

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

INTLIBRARY.

Definition at line 42 of file alloc.c.

Function Documentation

◆ base_copy()

Pbase base_copy ( Pbase  b)

Direct duplication.

The initial Pbase is assumed to be valid. Absolutely the same with base_dup, but base_up is the only function that maintains the old order. So recopy here for use with copy version including vect_copy, contrainte_copy, contraintes_copy, sc_copy (DN,24/6/02) Does not change the parameter. Did have a look at all copy version (DN,1/7/2002)

Definition at line 300 of file alloc.c.

301 {
302  Pbase n = BASE_NULLE, p = BASE_NULLE, r = BASE_NULLE, tmp = b;
303 
304  for (; tmp!=BASE_NULLE; tmp=tmp->succ)
305  {
306  n = (Pvecteur) MALLOC(sizeof(Svecteur),VECTEUR,"base_copy");
307  if (n == NULL) {
308  fprintf(stderr,"[base_copy] out of memory space\n");
309  abort();
310  }
311  var_of(n) = var_of(tmp);
312  val_of(n) = VALUE_ONE;
313  n->succ = NULL;
314  if (r==BASE_NULLE) r = n;
315  if (p!=BASE_NULLE) p->succ = n;
316  p = n;
317  }
318 
319  return r;
320 }
#define VALUE_ONE
#define abort()
Definition: misc-local.h:53
int fprintf()
test sc_min : ce test s'appelle par : programme fichier1.data fichier2.data ...
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)
struct Svecteur * Pvecteur
#define var_of(varval)
#define BASE_NULLE
MACROS SUR LES BASES.
#define VECTEUR
package sur les vecteurs creux et les bases
Definition: vecteur-local.h:51
#define MALLOC(s, t, f)
INTLIBRARY.
Definition: alloc.c:42

References abort, BASE_NULLE, fprintf(), MALLOC, Svecteur::succ, val_of, VALUE_ONE, var_of, and VECTEUR.

Referenced by base_append(), base_union(), build_integer_sc_nredund(), new_system_with_only_live_variable(), region_sc_projection_along_variables_ofl_ctrl(), regions_must_convex_hull(), sc_base_dup(), sc_copy(), sc_fm_project_variables(), sc_init_with_sc(), sc_of_constrs(), sc_projection_concat_proj_on_variables(), sc_safe_build_sc_nredund_1pass(), sc_safe_build_sc_nredund_2pass(), sc_safe_elim_redund(), sc_strong_normalize2(), sc_strong_normalize_and_check_feasibility2(), and sg_of_rays().

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

◆ base_dup()

Pbase base_dup ( Pbase  b)

Pbase base_dup(Pbase b) Note: this function changes the value of the pointer.

Use base_copy instead. Should become a link, not a function. For the moment, it's a function, because of the sc.h cannot be updated without installation, due to decision of integration of Janus or not? DN 12/5/03

return base_copy(b);

Definition at line 268 of file alloc.c.

269 {
270  /* return base_copy(b);*/
271 
272  Pbase n = BASE_NULLE, p = BASE_NULLE, r = BASE_NULLE;
273 
274  for (; b!=BASE_NULLE; b=b->succ)
275  {
276  n = (Pvecteur) MALLOC(sizeof(Svecteur),VECTEUR,"base_dup");
277  if (n == NULL) {
278  fprintf(stderr,"[base_dup] out of memory space\n");
279  abort();
280  }
281  var_of(n) = var_of(b);
282  val_of(n) = VALUE_ONE;
283  n->succ = NULL;
284  if (r==BASE_NULLE) r = n;
285  if (p!=BASE_NULLE) p->succ = n;
286  p = n;
287  }
288  return r;
289 }

References abort, BASE_NULLE, fprintf(), MALLOC, Svecteur::succ, val_of, VALUE_ONE, var_of, and VECTEUR.

Referenced by add_var_sup(), build_image_base(), build_sc_nredund_1pass_ofl_ctrl(), dependence_cone_positive(), dual_pivot(), loop_nest_to_wp65_code(), main(), movement_computation(), pip_solve(), pip_solve_min_with_big(), plint(), plreal(), primal(), primal_positive(), region_sc_minimal(), region_sc_projection_ofl_along_parameter(), regions_must_convex_hull(), sc_enveloppe_chernikova_ofl_ctrl(), sc_faisabilite_optim(), sc_projection_optim_along_vecteur_ofl(), sg_dup(), sort_psysteme(), syst_smith(), TestDependence(), transformer_convex_hulls(), transformer_derivative_fix_point(), transformer_equality_fix_point(), transformer_normalize(), transformer_projection_with_redundancy_elimination_and_check(), update_basis(), and var_ecart_sup().

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

◆ dbg_vect_rm()

void dbg_vect_rm ( Pvecteur  v,
char __attribute__((unused)) *  f 
)

void dbg_vect_rm(Pvecteur v, char * f): desallocation d'un vecteur avec marquage de la fonction provoquant la desallocation

Apparemment obsolete. RGSUSED

Definition at line 139 of file alloc.c.

141 {
142  Pvecteur v1,v2;
143  v1 = v;
144  while (v1!=NULL) {
145  v2 = v1->succ;
146  FREE((char *)v1,VECTEUR,f);
147  v1 = v2;
148  }
149 }
int f(int off1, int off2, int n, float r[n], float a[n], float b[n])
Definition: offsets.c:15
#define FREE(p, t, f)
Definition: alloc.c:43

References f(), FREE, Svecteur::succ, and VECTEUR.

Referenced by dbg_contrainte_rm(), dbg_ray_dte_rm(), and dbg_sommet_rm().

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

◆ vect_copy()

Pbase vect_copy ( Pvecteur  b)

direct duplication.

vect_dup() and vect_reversal() do the same thing : duplicate the vector with the reversal order. vect_copy duplicate the vector with the same order. in use of sc_copy. (DN,24/6/02) Does not change parameter b (DN,28/06/02)

Definition at line 240 of file alloc.c.

241 {
242  Pvecteur n = VECTEUR_NUL, p = VECTEUR_NUL, r = VECTEUR_NUL, tmp = b;
243 
244  for (; tmp!=VECTEUR_NUL; tmp=tmp->succ)
245  {
246  n = (Pvecteur) MALLOC(sizeof(Svecteur),VECTEUR,"vect_copy");
247  if (n == NULL) {
248  fprintf(stderr,"[vect_copy] out of memory space\n");
249  abort();
250  }
251  var_of(n) = var_of(tmp);
252  val_of(n) = val_of(tmp);
253  n->succ = NULL;
254  if (r==VECTEUR_NUL) r = n;
255  if (p!=VECTEUR_NUL) p->succ = n;
256  p = n;
257  }
258 
259  return r;
260 }
#define VECTEUR_NUL
DEFINITION DU VECTEUR NUL.

References abort, fprintf(), MALLOC, Svecteur::succ, val_of, var_of, VECTEUR, and VECTEUR_NUL.

Referenced by add_type_information(), c_convex_effects_on_actual_parameter_forward_translation(), check_positive_dependence(), contrainte_copy(), contrainte_parallele(), find_motif(), find_pattern(), main(), make_bound_expression(), print_call_precondition(), print_loopnest_dependence_cone(), sc_copy(), sc_elim_double_constraints(), sc_minmax_of_variable(), sc_projection_concat_proj_on_variables(), small_positive_slope_reduce_coefficients_with_bounding_box(), vect_gen_copy_tree(), xml_GlobalVariables(), and xml_LocalVariables().

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

◆ vect_dup()

Pvecteur vect_dup ( Pvecteur  v_in)

Pvecteur vect_dup(Pvecteur v_in): duplication du vecteur v_in; allocation de et copie dans v_out;.

end of vecteur-local.h

allocate v_out; v_out := v_in;

Parameters
v_in_in

Definition at line 51 of file alloc.c.

53 {
54  Pvecteur v_out;
55  Pvecteur v;
56 
57  v_out = NULL;
58  for(v=v_in; v!=NULL; v=v->succ) {
59  v_out = vect_chain(v_out,var_of(v),val_of(v));
60  }
61 
62  return v_out;
63 }
Pvecteur vect_chain(Pvecteur v_in, Variable var, Value coeff)
package vecteur routines internes au package
Definition: private.c:69

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

Referenced by __attribute__(), add_affine_bound_conditions(), adg_sc_dup(), affine_to_transformer(), array_indices_communication(), block_to_complexity(), bounds_equal_p(), broadcast_of_dataflow(), build_list_of_min(), build_sc_machine(), build_sc_with_several_uniform_ref(), build_third_comb(), calculate_delay(), change_base_in_sc(), check_range_wrt_precondition(), complex_bound_computation(), constraint_to_bound(), contrainte_dup(), CopyAccVec(), create_farkas_poly(), dependence_system_add_lci_and_di(), dj_simple_inegs_to_eg(), dj_system_complement(), dup_list_of_Pvecteur(), dup_vv(), elim_var_with_eg(), eq_in_ineq(), expression_to_affine(), final_statement_to_complexity_evaluation(), find_eg(), find_pattern(), formal_and_actual_parameters_association(), free_guards(), fusion_buffer(), gcd_and_constant_dependence_test(), get_bounds_expression(), gomory_trait_eq(), include_trans_on_LC_in_ref(), integer_divide_to_transformer(), integer_power_to_transformer(), integer_right_shift_to_transformer(), loop_executed_approximation(), loop_index_domaine_to_contrainte(), loop_nest_to_wp65_code(), lower_bound_generation(), make_constraint_expression(), make_context_of_loop(), make_datum_movement(), make_dual(), make_load_blocks(), make_movements_loop_body_wp65(), make_reindex(), make_store_blocks(), matrices_to_loop_sc(), matrices_to_sc(), matrix_to_system(), MaxBoundary(), MinBoundary(), monome_dup(), monome_monome_div(), monome_monome_mult(), movement_computation(), my_vect_substract(), negate_expression(), normalize_intrinsic(), pa_path_to_few_disjunct_ofl_ctrl(), parallel_tiling(), Pcontrainte_separate_on_vars(), pivoter_pas(), plc_elim_var_with_eg(), plint_degen(), Pvecteur_to_assign_statement(), ray_dte_dup(), sc_dup1(), sc_find_equalities(), sc_image_computation(), sc_lexicographic_sort(), sc_projection_optim_along_vecteur_ofl(), sc_supress_same_constraints(), sc_to_tableau(), separate_variables(), set_dimensions_of_local_variable_family(), simple_affine_to_transformer(), simplify_deducable_variables(), simplify_dimension(), simplify_minmax_contrainte(), simplify_predicate(), sommet_dup(), substitute_and_create(), tile_hyperplane_constraints(), tile_membership_constraints(), Tiling_buffer_allocation(), tiling_transformation(), transformer_add_integer_relation_information(), transformer_add_loop_index_initialization(), transformer_equalities_add(), transformer_pattern_fix_point(), translate_global_values(), update_basis(), upper_bound_generation(), vect_add(), vect_add_first(), vect_change_base(), vect_cl_ofl_ctrl(), vect_del_var(), vect_extract(), vect_sort(), vect_substract(), xml_Compute_and_Need(), xml_Pattern_Paving(), and xml_Region_Range().

+ Here is the call graph for this function:

◆ vect_make()

Pvecteur vect_make ( Pvecteur  v,
Variable  var,
Value  val,
  ... 
)

Pvecteur vect_make(v, [var, val,]* 0, val) Pvecteur v; // may be NULL, use assigne anyway Variable var; Value val;.

Builds a vector from the list of arguments, by successive additions. ends when a 0 Variable (that is TCST!) is encountered.

Because of the var val order, this function cannot be called directly with a va_list, but (va_list, 0) should be used, since the val argument is expected, read and used anyway.

CAUTION: the initial vector is modified by the process!

Parameters
varar
valal

Definition at line 165 of file alloc.c.

166 {
167  va_list the_args;
168 
169  // handle fist argument
170  vect_add_elem(&v, var, val);
171 
172  // get others
173  va_start(the_args, val);
174 
175  while (var != (Variable) 0)
176  {
177  var = va_arg(the_args, Variable);
178  val = va_arg(the_args, Value);
179  vect_add_elem(&v, var, val);
180  }
181 
182  va_end (the_args);
183  return v;
184 }
int Value
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
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 vect_add_elem().

Referenced by add_type_information(), align_check(), buffer_acces(), call_rwt(), compute_entity_to_declaration_constraints(), compute_x_and_y_bounds(), convex_in_effect_loop_range_fix(), fortran_data_to_prec_for_variables(), fortran_user_function_call_to_transformer(), fusion(), fusion_buffer(), generate_system_for_equal_variables(), generate_work_sharing_system(), Hierarchical_tiling(), hpfc_compute_align_constraints(), hpfc_compute_distribute_constraints(), hpfc_compute_entity_to_new_declaration(), hpfc_compute_unicity_constraints(), interlaced_basic_workchunk_regions_p(), loop_basic_workchunk_to_workchunk(), make_loop_indice_equation(), processor_loop(), sc_image_computation(), sc_minmax_of_variables(), small_positive_slope_reduce_coefficients_with_bounding_box(), Tiling2_buffer(), Tiling_buffer_allocation(), transformer_derivative_constraints(), transformer_derivative_fix_point(), and transformer_list_generic_transitive_closure().

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

◆ vect_make_1D()

Pvecteur vect_make_1D ( Value  a,
Variable  x,
Value  b 
)

Generate a sparse vector a x + b TCST.

Definition at line 226 of file alloc.c.

227 {
228  Pvecteur v = vect_new(x, a);
229  vect_add_elem(&v, TCST, b);
230  return v;
231 }
static char * x
Definition: split_file.c:159
#define TCST
VARIABLE REPRESENTANT LE TERME CONSTANT.
Pvecteur vect_new(Variable var, Value coeff)
Pvecteur vect_new(Variable var,Value coeff): allocation d'un vecteur colineaire au vecteur de base va...
Definition: alloc.c:110

References TCST, vect_add_elem(), vect_new(), and x.

Referenced by contrainte_make_1D().

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

◆ vect_make_dense()

Pvecteur vect_make_dense ( Pbase  b,
Value  val,
  ... 
)

Allocate a new vector v whose coefficient are given by the list of values ad whose dimension is given by b.

The number of constant values passed as argument is supposed to be equal to the dimension of b.

Note: 0 is a normal value. I see no way to mark the last argument.

FI: I add this function to check under gdb that a given point belongs to a constraint system. The manual verification is tedious and error prone. This is done for debugging in Linear/C3 Library the linked_regions bug.

handle first argument - the first element of a basis has rank 1

get others

Parameters
valal

Definition at line 198 of file alloc.c.

199 {
200  va_list the_args;
201  Pvecteur v = VECTEUR_NUL;
202  Variable var;
203  int dim = (int) base_dimension(b);
204  int i = 0; // current dimension
205 
206  /* handle first argument - the first element of a basis has rank 1 */
207  var = variable_of_rank(b,1);
208  vect_add_elem(&v, var, val);
209 
210  /* get others */
211  va_start (the_args, val);
212 
213  for(i=2;i<=dim;i++) // i <= dim because of the way dimensions are counted
214  {
215  var = variable_of_rank(b,i);
216  val = va_arg(the_args, Value);
217  vect_add_elem(&v, var, val);
218  }
219 
220  va_end (the_args);
221 
222  return v;
223 }
void const char const char const int
Variable variable_of_rank()
#define base_dimension(b)

References base_dimension, int, variable_of_rank(), vect_add_elem(), and VECTEUR_NUL.

Referenced by build_sc_nredund_2pass_ofl_ctrl().

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

◆ vect_new()

Pvecteur vect_new ( Variable  var,
Value  coeff 
)

Pvecteur vect_new(Variable var,Value coeff): allocation d'un vecteur colineaire au vecteur de base var et de coefficient coeff (i.e.

creation d'un nouveau vecteur ne comportant qu'un seul couple (var,coeff))

  -->

coeff var

Pourrait etre remplace par un vect_chain(NULL,,)

Modifications:

  • a 0 coeff generates a null vector; Francois Irigoin, 26 March 1991

fprintf(stderr, "%10.3f MB", (sbrk(0) - etext)/(double)(1 << 20)); // not portable

xit(-1);

Parameters
varar
coeffoeff

Definition at line 110 of file alloc.c.

111 {
112  Pvecteur v;
113 
114  if(coeff!=0) {
115  v = (Pvecteur) MALLOC(sizeof(Svecteur),VECTEUR,"vect_new");
116  if (v == NULL) {
117  (void) fprintf(stderr,"vect_new: Out of memory space\n");
118  /* fprintf(stderr, "%10.3f MB",
119  (sbrk(0) - etext)/(double)(1 << 20)); // not portable */
120  abort();
121  /*exit(-1);*/
122  }
123  var_of(v) = var;
124  val_of(v) = coeff;
125  v->succ = NULL;
126  }
127  else
128  v = NULL;
129 
130  return v;
131 }

References abort, fprintf(), MALLOC, Svecteur::succ, val_of, var_of, and VECTEUR.

Referenced by add_elem_to_list_of_Pvecteur(), add_equivalence_equality(), add_loop_index_exit_value(), add_type_information(), add_x_list(), adg_dataflowgraph(), adg_dataflowgraph_with_extremities(), adg_get_predicate_of_loops(), adg_max_of_leaves(), affine_to_transformer(), ajout_dte(), align_check(), apply_farkas(), array_indices_communication(), array_scalar_access_to_bank_communication(), base_add_variable(), base_complete(), binary_to_normalized(), bitwise_xor_to_transformer(), broadcast_of_dataflow(), build_esv_list(), build_image_base(), build_sc_machine(), build_third_comb(), c_convex_effects_on_actual_parameter_forward_translation(), calculate_delay(), comp_exec_domain(), completer_base(), config_vecteur(), converti_psysmin_psysmax(), cout_nul(), create_farkas_poly(), dependence_cone_positive(), dj_system_complement(), do_solve_hardware_constraints_on_nb_proc(), do_terapix_warmup_patching(), elim_var_with_eg(), entity_list_to_base(), expression_equal_in_context_p(), expression_flt(), expression_less_than_in_context(), expression_multiply_sizeof_to_transformer(), find_implicit_equation(), find_vbase(), fonct_max(), fonct_max_all(), fonct_max_d(), fonct_min(), fonct_min_all(), fonct_min_d(), fonct_read(), generate_one_message(), generate_work_sharing_system(), generic_abs_to_transformer(), generic_equality_to_transformer(), generic_minmax_to_transformer(), hpfc_compute_distribute_constraints(), hpfc_compute_entity_to_new_declaration(), iabs_to_transformer(), include_parameters_in_sc(), integer_left_shift_to_transformer(), integer_minmax_to_transformer(), integer_multiply_to_transformer(), integer_power_to_transformer(), local_tile_constraints(), logical_binary_function_to_transformer(), logical_binary_operation_to_transformer(), logical_constant_to_transformer(), logical_unary_operation_to_transformer(), loop_basic_workchunk_to_workchunk(), loop_bound_evaluation_to_transformer(), loop_bounds_to_tile_bounds(), loop_index_domaine_to_contrainte(), loop_regions_normalize(), lvbase_add(), make_datum_movement(), make_dual(), make_load_blocks(), make_monome(), make_movement_scalar_wp65(), make_movements_loop_body_wp65(), make_store_blocks(), mat_sys_conv(), matrices_to_constraints(), matrices_to_constraints_with_sym_cst(), matrices_to_contraintes_with_sym_cst(), matrices_to_loop_sc(), matrices_to_sc(), matrix_to_system(), mk_rn(), monome_del_var(), monome_monome_div(), monome_monome_mult(), movement_computation(), my_matrices_to_constraints_with_sym_cst(), my_matrices_to_constraints_with_sym_cst_2(), my_substitute_var_with_vec(), normalize_constant(), normalize_reference(), NormalizeCall(), NormalizeConstant(), NormalizeReference(), old_prototype_factorize(), one_message_guards_and_neighbour(), one_receive_message(), pa_path_to_few_disjunct_ofl_ctrl(), plc_elim_var_with_eg(), plc_make_vvs_with_vector(), plint_degen(), polynome_constant_p(), polynome_roots(), polynome_TCST(), prepare_reindexing(), prototype_factorize(), pu_matrices_to_contraintes(), put_source_ind(), re_do_it(), reference_conversion_computation(), relation_to_transformer(), remove_variables_if_possible(), sc_base_add_variable(), sc_empty(), sc_find_equalities(), sc_image_computation(), sc_multiply_constant_terms(), sc_to_vvs(), simple_addition_to_transformer(), simple_affine_to_transformer(), simplify_float_constraint(), sort_tile_indices(), substitute_var_with_vec(), top_down_abc_dimension(), 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_sign_information(), transformer_convex_hulls(), transformer_logical_inequalities_add(), translate_global_value(), valuer(), var_pivotd(), var_pivots(), var_posit(), variables_in_declaration_list(), vect_chain(), vect_make_1D(), vect_make_line(), vect_read(), vecteur_mult(), vecteur_of_zvec(), xml_Chain_Graph(), xml_LocalVariables(), and xml_Pattern_Paving().

+ Here is the call graph for this function:

◆ vect_rm()

void vect_rm ( Pvecteur  v)

void vect_rm(Pvecteur v): desallocation des couples de v;

Attention! La procedure appelante doit penser a modifier la valeur de v apres l'appel pour ne pas pointer dans le vide: vect_rm(v); v = NULL;

Il vaudrait mieux que vect_rm retourne NULL et qu'on puisse ecrire: v = vect_rm(v); ou que vect_rm prenne un Pvecteur * en argument: vect_rm(&v);

Definition at line 78 of file alloc.c.

80 {
81  Pvecteur elem = NULL;
82  Pvecteur nelem = NULL;
83 
84  for(elem = v; elem!=NULL; elem = nelem) {
85  nelem = elem->succ;
86  free((char *)elem);
87  }
88  /*
89  while (v != NULL) {
90  Pvecteur nv = v->succ;
91  free((char *) v);
92  v = nv;
93  }
94  */
95 }
void free(void *)

References free(), and Svecteur::succ.

Referenced by add_affine_bound_conditions(), add_reference_information(), add_type_information(), affine_expression_of_loop_index_p(), affine_to_transformer(), append_eg(), bounds_equal_p(), build_and_test_dependence_context(), build_convex_constraints_from_vertices(), build_image_base(), c_convex_effects_on_actual_parameter_forward_translation(), cell_reference_sc_exact_projection_along_variable(), check_positive_dependence(), constraints_keep_invariants_only(), constraints_to_loop_bound(), contrainte_free(), contrainte_parallele(), contrainte_reversal(), contrainte_subst_ofl_ctrl(), cout_nul(), dependence_system_add_lci_and_di(), dj_system_complement(), do_group_statement_constant(), do_solve_hardware_constraints_on_volume(), do_terapix_warmup_patching(), dual(), dual_positive(), elim_var_with_eg(), eq_in_ineq(), eq_set_vect_nul(), erase_trivial_ineg(), expression_equal_in_context_p(), expression_flt(), expression_less_than_in_context(), find_eg(), find_motif(), find_pattern(), find_vbase(), free_vector_list(), FreeNormalized(), gcd_and_constant_dependence_test(), hpfc_integer_constant_expression_p(), invariant_vector_p(), is_inferior_monome(), lvbase_ote_no_ligne(), make_bound_expression(), make_vecteur_expression(), monome_monome_mult(), monome_rm(), my_substitute_var_with_vec(), my_system_remove_variables(), outliner_smart_references_computation(), pa_path_to_few_disjunct_ofl_ctrl(), Pcontrainte_to_expression_list(), plc_elim_var_with_eg(), plint(), plint_degen(), plreal(), polynome_constant_p(), polynome_TCST(), polynome_used_var(), primal(), primal_positive(), print_call_precondition(), print_loopnest_dependence_cone(), Pvecteur_to_assign_statement(), ray_dte_rm(), region_exact_projection_along_parameters(), region_exact_projection_along_variable(), remove_variables_if_possible(), sc_bounded_normalization(), sc_concatenate(), sc_elim_db_constraints(), sc_elim_double_constraints(), sc_elim_redund_with_first_ofl_ctrl(), sc_find_equalities(), sc_free1(), sc_kill_db_eg(), sc_lexicographic_sort(), sc_minmax_of_variable(), sc_minmax_of_variable_optim(), sc_proj_optim_on_di_ofl(), sc_projection_ofl_along_list_of_variables(), sc_projection_optim_along_vecteur_ofl(), sc_rm(), sc_safe_elim_db_constraints(), sc_safe_kill_db_eg(), sc_simplex_feasibility_ofl_ctrl_fixprec(), sc_substitute_dimension(), sc_supress_parallel_redund_constraints(), sc_supress_same_constraints(), sc_transform_ineg_in_eg(), shift_expression_of_loop_index_p(), simple_affine_to_transformer(), simplify_constraint_with_bounding_box(), simplify_float_constraint(), simplify_float_constraint_system(), sl_fprint_tab(), sommet_rm(), sommets_rm(), substitute_and_create(), substitute_var_with_vec(), top_down_abc_dimension(), transformer_add_integer_relation_information(), transformer_projection_with_redundancy_elimination_and_check(), translate_to_module_frame(), try_reorder_expression_call(), uniform_dependence_p(), valuer(), var_pivotd(), var_pivots(), vect_gen_free(), vect_multiply(), vect_product(), vect_subst(), xml_Chain_Graph(), xml_Compute_and_Need(), xml_GlobalVariables(), and xml_LocalVariables().

+ Here is the call graph for this function: