MLIR
20.0.0git
|
An IntegerRelation represents the set of points from a PresburgerSpace that satisfy a list of affine constraints. More...
#include "mlir/Analysis/Presburger/IntegerRelation.h"
Classes | |
struct | CountsSnapshot |
The struct CountsSnapshot stores the count of each VarKind, and also of each constraint type. More... | |
Public Types | |
enum class | Kind { IntegerRelation , IntegerPolyhedron , FlatLinearConstraints , FlatLinearValueConstraints , FlatAffineValueConstraints , FlatAffineRelation } |
All derived classes of IntegerRelation. More... | |
Public Member Functions | |
IntegerRelation (unsigned numReservedInequalities, unsigned numReservedEqualities, unsigned numReservedCols, const PresburgerSpace &space) | |
Constructs a relation reserving memory for the specified number of constraints and variables. More... | |
IntegerRelation (const PresburgerSpace &space) | |
Constructs a relation with the specified number of dimensions and symbols. More... | |
virtual | ~IntegerRelation ()=default |
virtual Kind | getKind () const |
Return the kind of this IntegerRelation. More... | |
std::unique_ptr< IntegerRelation > | clone () const |
const PresburgerSpace & | getSpace () const |
Returns a reference to the underlying space. More... | |
void | setSpace (const PresburgerSpace &oSpace) |
Set the space to oSpace , which should have the same number of ids as the current space. More... | |
void | setSpaceExceptLocals (const PresburgerSpace &oSpace) |
Set the space to oSpace , which should not have any local ids. More... | |
void | setId (VarKind kind, unsigned i, Identifier id) |
Set the identifier for the ith variable of the specified kind of the IntegerRelation's PresburgerSpace. More... | |
void | resetIds () |
ArrayRef< Identifier > | getIds (VarKind kind) |
Get the identifiers for the variables of specified varKind. More... | |
PresburgerSpace | getSpaceWithoutLocals () const |
Returns a copy of the space without locals. More... | |
void | append (const IntegerRelation &other) |
Appends constraints from other into this . More... | |
IntegerRelation | intersect (IntegerRelation other) const |
Return the intersection of the two relations. More... | |
bool | isEqual (const IntegerRelation &other) const |
Return whether this and other are equal. More... | |
bool | isObviouslyEqual (const IntegerRelation &other) const |
Perform a quick equality check on this and other . More... | |
bool | isSubsetOf (const IntegerRelation &other) const |
Return whether this is a subset of the given IntegerRelation. More... | |
DynamicAPInt | atEq (unsigned i, unsigned j) const |
Returns the value at the specified equality row and column. More... | |
int64_t | atEq64 (unsigned i, unsigned j) const |
The same, but casts to int64_t. More... | |
DynamicAPInt & | atEq (unsigned i, unsigned j) |
DynamicAPInt | atIneq (unsigned i, unsigned j) const |
Returns the value at the specified inequality row and column. More... | |
int64_t | atIneq64 (unsigned i, unsigned j) const |
The same, but casts to int64_t. More... | |
DynamicAPInt & | atIneq (unsigned i, unsigned j) |
unsigned | getNumConstraints () const |
unsigned | getNumDomainVars () const |
unsigned | getNumRangeVars () const |
unsigned | getNumSymbolVars () const |
unsigned | getNumLocalVars () const |
unsigned | getNumDimVars () const |
unsigned | getNumDimAndSymbolVars () const |
unsigned | getNumVars () const |
unsigned | getNumCols () const |
Returns the number of columns in the constraint system. More... | |
unsigned | getNumEqualities () const |
unsigned | getNumInequalities () const |
unsigned | getNumReservedEqualities () const |
unsigned | getNumReservedInequalities () const |
ArrayRef< DynamicAPInt > | getEquality (unsigned idx) const |
ArrayRef< DynamicAPInt > | getInequality (unsigned idx) const |
SmallVector< int64_t, 8 > | getEquality64 (unsigned idx) const |
The same, but casts to int64_t. More... | |
SmallVector< int64_t, 8 > | getInequality64 (unsigned idx) const |
IntMatrix | getInequalities () const |
unsigned | getNumVarKind (VarKind kind) const |
Get the number of vars of the specified kind. More... | |
unsigned | getVarKindOffset (VarKind kind) const |
Return the index at which the specified kind of vars starts. More... | |
unsigned | getVarKindEnd (VarKind kind) const |
Return the index at Which the specified kind of vars ends. More... | |
unsigned | getVarKindOverlap (VarKind kind, unsigned varStart, unsigned varLimit) const |
Get the number of elements of the specified kind in the range [varStart, varLimit). More... | |
VarKind | getVarKindAt (unsigned pos) const |
Return the VarKind of the var at the specified position. More... | |
CountsSnapshot | getCounts () const |
void | truncate (const CountsSnapshot &counts) |
virtual unsigned | insertVar (VarKind kind, unsigned pos, unsigned num=1) |
Insert num variables of the specified kind at position pos . More... | |
unsigned | appendVar (VarKind kind, unsigned num=1) |
Append num variables of the specified kind after the last variable of that kind. More... | |
void | addInequality (ArrayRef< DynamicAPInt > inEq) |
Adds an inequality (>= 0) from the coefficients specified in inEq . More... | |
void | addInequality (ArrayRef< int64_t > inEq) |
void | addEquality (ArrayRef< DynamicAPInt > eq) |
Adds an equality from the coefficients specified in eq . More... | |
void | addEquality (ArrayRef< int64_t > eq) |
virtual void | eliminateRedundantLocalVar (unsigned posA, unsigned posB) |
Eliminate the posB^th local variable, replacing every instance of it with the posA^th local variable. More... | |
void | removeVar (VarKind kind, unsigned pos) |
Removes variables of the specified kind with the specified pos (or within the specified range) from the system. More... | |
virtual void | removeVarRange (VarKind kind, unsigned varStart, unsigned varLimit) |
void | removeVar (unsigned pos) |
Removes the specified variable from the system. More... | |
void | removeEquality (unsigned pos) |
void | removeInequality (unsigned pos) |
void | removeEqualityRange (unsigned start, unsigned end) |
Remove the (in)equalities at positions [start, end). More... | |
void | removeInequalityRange (unsigned start, unsigned end) |
MaybeOptimum< SmallVector< Fraction, 8 > > | findRationalLexMin () const |
Get the lexicographically minimum rational point satisfying the constraints. More... | |
MaybeOptimum< SmallVector< DynamicAPInt, 8 > > | findIntegerLexMin () const |
Same as above, but returns lexicographically minimal integer point. More... | |
virtual void | swapVar (unsigned posA, unsigned posB) |
Swap the posA^th variable with the posB^th variable. More... | |
void | clearConstraints () |
Removes all equalities and inequalities. More... | |
void | setAndEliminate (unsigned pos, ArrayRef< DynamicAPInt > values) |
Sets the values.size() variables starting at po s to the specified values and removes them. More... | |
void | setAndEliminate (unsigned pos, ArrayRef< int64_t > values) |
virtual void | clearAndCopyFrom (const IntegerRelation &other) |
Replaces the contents of this IntegerRelation with other . More... | |
void | getLowerAndUpperBoundIndices (unsigned pos, SmallVectorImpl< unsigned > *lbIndices, SmallVectorImpl< unsigned > *ubIndices, SmallVectorImpl< unsigned > *eqIndices=nullptr, unsigned offset=0, unsigned num=0) const |
Gather positions of all lower and upper bounds of the variable at pos , and optionally any equalities on it. More... | |
bool | isEmpty () const |
Checks for emptiness by performing variable elimination on all variables, running the GCD test on each equality constraint, and checking for invalid constraints. More... | |
bool | isObviouslyEmpty () const |
Performs GCD checks and invalid constraint checks. More... | |
bool | isEmptyByGCDTest () const |
Runs the GCD test on all equality constraints. More... | |
bool | isIntegerEmpty () const |
Returns true if the set of constraints is found to have no solution, false if a solution exists. More... | |
IntMatrix | getBoundedDirections () const |
Returns a matrix where each row is a vector along which the polytope is bounded. More... | |
std::optional< SmallVector< DynamicAPInt, 8 > > | findIntegerSample () const |
Find an integer sample point satisfying the constraints using a branch and bound algorithm with generalized basis reduction, with some additional processing using Simplex for unbounded sets. More... | |
std::optional< DynamicAPInt > | computeVolume () const |
Compute an overapproximation of the number of integer points in the relation. More... | |
bool | containsPoint (ArrayRef< DynamicAPInt > point) const |
Returns true if the given point satisfies the constraints, or false otherwise. More... | |
bool | containsPoint (ArrayRef< int64_t > point) const |
std::optional< SmallVector< DynamicAPInt, 8 > > | containsPointNoLocal (ArrayRef< DynamicAPInt > point) const |
Given the values of non-local vars, return a satisfying assignment to the local if one exists, or an empty optional otherwise. More... | |
std::optional< SmallVector< DynamicAPInt, 8 > > | containsPointNoLocal (ArrayRef< int64_t > point) const |
DivisionRepr | getLocalReprs (std::vector< MaybeLocalRepr > *repr=nullptr) const |
Returns a DivisonRepr representing the division representation of local variables in the constraint system. More... | |
void | addBound (BoundType type, unsigned pos, const DynamicAPInt &value) |
Adds a constant bound for the specified variable. More... | |
void | addBound (BoundType type, unsigned pos, int64_t value) |
void | addBound (BoundType type, ArrayRef< DynamicAPInt > expr, const DynamicAPInt &value) |
Adds a constant bound for the specified expression. More... | |
void | addBound (BoundType type, ArrayRef< int64_t > expr, int64_t value) |
void | addLocalFloorDiv (ArrayRef< DynamicAPInt > dividend, const DynamicAPInt &divisor) |
Adds a new local variable as the floordiv of an affine function of other variables, the coefficients of which are provided in dividend and with respect to a positive constant divisor . More... | |
void | addLocalFloorDiv (ArrayRef< int64_t > dividend, int64_t divisor) |
void | projectOut (unsigned pos, unsigned num) |
Projects out (aka eliminates) num variables starting at position pos . More... | |
void | projectOut (unsigned pos) |
LogicalResult | constantFoldVar (unsigned pos) |
Tries to fold the specified variable to a constant using a trivial equality detection; if successful, the constant is substituted for the variable everywhere in the constraint system and then removed from the system. More... | |
void | constantFoldVarRange (unsigned pos, unsigned num) |
This method calls constantFoldVar for the specified range of variables, num variables starting at position pos . More... | |
LogicalResult | unionBoundingBox (const IntegerRelation &other) |
Updates the constraints to be the smallest bounding (enclosing) box that contains the points of this set and that of other , with the symbols being treated specially. More... | |
std::optional< DynamicAPInt > | getConstantBoundOnDimSize (unsigned pos, SmallVectorImpl< DynamicAPInt > *lb=nullptr, DynamicAPInt *boundFloorDivisor=nullptr, SmallVectorImpl< DynamicAPInt > *ub=nullptr, unsigned *minLbPos=nullptr, unsigned *minUbPos=nullptr) const |
Returns the smallest known constant bound for the extent of the specified variable (pos^th), i.e., the smallest known constant that is greater than or equal to 'exclusive upper bound' - 'lower bound' of the variable. More... | |
std::optional< int64_t > | getConstantBoundOnDimSize64 (unsigned pos, SmallVectorImpl< int64_t > *lb=nullptr, int64_t *boundFloorDivisor=nullptr, SmallVectorImpl< int64_t > *ub=nullptr, unsigned *minLbPos=nullptr, unsigned *minUbPos=nullptr) const |
The same, but casts to int64_t. More... | |
std::optional< DynamicAPInt > | getConstantBound (BoundType type, unsigned pos) const |
Returns the constant bound for the pos^th variable if there is one; std::nullopt otherwise. More... | |
std::optional< int64_t > | getConstantBound64 (BoundType type, unsigned pos) const |
The same, but casts to int64_t. More... | |
void | removeIndependentConstraints (unsigned pos, unsigned num) |
Removes constraints that are independent of (i.e., do not have a coefficient) variables in the range [pos, pos + num). More... | |
bool | isHyperRectangular (unsigned pos, unsigned num) const |
Returns true if the set can be trivially detected as being hyper-rectangular on the specified contiguous set of variables. More... | |
void | removeTrivialRedundancy () |
Removes duplicate constraints, trivially true constraints, and constraints that can be detected as redundant as a result of differing only in their constant term part. More... | |
void | removeRedundantInequalities () |
A more expensive check than removeTrivialRedundancy to detect redundant inequalities. More... | |
void | removeRedundantConstraints () |
Removes redundant constraints using Simplex. More... | |
void | removeDuplicateDivs () |
void | simplify () |
Simplify the constraint system by removing canonicalizing constraints and removing redundant constraints. More... | |
void | convertVarKind (VarKind srcKind, unsigned varStart, unsigned varLimit, VarKind dstKind, unsigned pos) |
Converts variables of kind srcKind in the range [varStart, varLimit) to variables of kind dstKind. More... | |
void | convertVarKind (VarKind srcKind, unsigned varStart, unsigned varLimit, VarKind dstKind) |
void | convertToLocal (VarKind kind, unsigned varStart, unsigned varLimit) |
void | mergeAndAlignSymbols (IntegerRelation &other) |
Merge and align symbol variables of this and other with respect to identifiers. More... | |
unsigned | mergeLocalVars (IntegerRelation &other) |
Adds additional local vars to the sets such that they both have the union of the local vars in each set, without changing the set of points that lie in this and other . More... | |
bool | hasOnlyDivLocals () const |
Check whether all local ids have a division representation. More... | |
void | setDimSymbolSeparation (unsigned newSymbolCount) |
Changes the partition between dimensions and symbols. More... | |
IntegerPolyhedron | getDomainSet () const |
Return a set corresponding to all points in the domain of the relation. More... | |
IntegerPolyhedron | getRangeSet () const |
Return a set corresponding to all points in the range of the relation. More... | |
void | intersectDomain (const IntegerPolyhedron &poly) |
Intersect the given poly with the domain in-place. More... | |
void | intersectRange (const IntegerPolyhedron &poly) |
Intersect the given poly with the range in-place. More... | |
void | inverse () |
Invert the relation i.e., swap its domain and range. More... | |
void | compose (const IntegerRelation &rel) |
Let the relation this be R1, and the relation rel be R2. More... | |
void | applyDomain (const IntegerRelation &rel) |
Given a relation rel , apply the relation to the domain of this relation. More... | |
void | applyRange (const IntegerRelation &rel) |
Given a relation rel , apply the relation to the range of this relation. More... | |
void | mergeAndCompose (const IntegerRelation &other) |
Given a relation other: (A -> B) , this operation merges the symbol and local variables and then takes the composition of other on this: (B -> C) . More... | |
PresburgerRelation | computeReprWithOnlyDivLocals () const |
Compute an equivalent representation of the same set, such that all local vars in all disjuncts have division representations. More... | |
SymbolicLexOpt | findSymbolicIntegerLexMin () const |
Compute the symbolic integer lexmin of the relation. More... | |
SymbolicLexOpt | findSymbolicIntegerLexMax () const |
Same as findSymbolicIntegerLexMin but produces lexmax instead of lexmin. More... | |
PresburgerRelation | subtract (const PresburgerRelation &set) const |
Return the set difference of this set and the given set, i.e., return this \ set . More... | |
void | removeTrivialEqualities () |
bool | isFullDim () |
void | print (raw_ostream &os) const |
void | dump () const |
Static Public Member Functions | |
static IntegerRelation | getUniverse (const PresburgerSpace &space) |
Return a system with no constraints, i.e., one which is satisfied by all points. More... | |
static IntegerRelation | getEmpty (const PresburgerSpace &space) |
Return an empty system containing an invalid equation 0 = 1. More... | |
static bool | classof (const IntegerRelation *cst) |
Protected Member Functions | |
bool | hasInvalidConstraint () const |
Checks all rows of equality/inequality constraints for trivial contradictions (for example: 1 == 0, 0 >= 1), which may have surfaced after elimination. More... | |
template<bool isLower> | |
std::optional< DynamicAPInt > | computeConstantLowerOrUpperBound (unsigned pos) |
Returns the constant lower bound if isLower is true, and the upper bound if isLower is false. More... | |
template<bool isLower> | |
std::optional< int64_t > | computeConstantLowerOrUpperBound64 (unsigned pos) |
The same, but casts to int64_t. More... | |
LogicalResult | gaussianEliminateVar (unsigned position) |
Eliminates a single variable at position from equality and inequality constraints. More... | |
void | removeRedundantLocalVars () |
Removes local variables using equalities. More... | |
unsigned | gaussianEliminateVars (unsigned posStart, unsigned posLimit) |
Eliminates variables from equality and inequality constraints in column range [posStart, posLimit). More... | |
bool | gaussianEliminate () |
Perform a Gaussian elimination operation to reduce all equations to standard form. More... | |
virtual void | fourierMotzkinEliminate (unsigned pos, bool darkShadow=false, bool *isResultIntegerExact=nullptr) |
Eliminates the variable at the specified position using Fourier-Motzkin variable elimination, but uses Gaussian elimination if there is an equality involving that variable. More... | |
void | gcdTightenInequalities () |
Tightens inequalities given that we are dealing with integer spaces. More... | |
void | normalizeConstraintsByGCD () |
Normalized each constraints by the GCD of its coefficients. More... | |
bool | findConstraintWithNonZeroAt (unsigned colIdx, bool isEq, unsigned *rowIdx) const |
Searches for a constraint with a non-zero coefficient at colIdx in equality (isEq=true) or inequality (isEq=false) constraints. More... | |
bool | isColZero (unsigned pos) const |
Returns true if the pos^th column is all zero for both inequalities and equalities. More... | |
bool | removeDuplicateConstraints () |
Checks for identical inequalities and eliminates redundant inequalities. More... | |
virtual bool | hasConsistentState () const |
Returns false if the fields corresponding to various variable counts, or equality/inequality buffer sizes aren't consistent; true otherwise. More... | |
virtual void | printSpace (raw_ostream &os) const |
Prints the number of constraints, dimensions, symbols and locals in the IntegerRelation. More... | |
void | removeVarRange (unsigned varStart, unsigned varLimit) |
Removes variables in the column range [varStart, varLimit), and copies any remaining valid data into place, updates member variables, and resizes arrays as needed. More... | |
void | truncateVarKind (VarKind kind, unsigned num) |
Truncate the vars of the specified kind to the specified number by dropping some vars at the end. More... | |
void | truncateVarKind (VarKind kind, const CountsSnapshot &counts) |
Truncate the vars to the number in the space of the specified CountsSnapshot. More... | |
Protected Attributes | |
PresburgerSpace | space |
IntMatrix | equalities |
Coefficients of affine equalities (in == 0 form). More... | |
IntMatrix | inequalities |
Coefficients of affine inequalities (in >= 0 form). More... | |
Static Protected Attributes | |
constexpr static unsigned | kExplosionFactor = 32 |
A parameter that controls detection of an unrealistic number of constraints. More... | |
An IntegerRelation represents the set of points from a PresburgerSpace that satisfy a list of affine constraints.
Affine constraints can be inequalities or equalities in the form:
Inequality: c_0*x_0 + c_1*x_1 + .... + c_{n-1}*x_{n-1} + c_n >= 0 Equality : c_0*x_0 + c_1*x_1 + .... + c_{n-1}*x_{n-1} + c_n == 0
where c_0, c_1, ..., c_n are integers and n is the total number of variables in the space.
Such a relation corresponds to the set of integer points lying in a convex polyhedron. For example, consider the relation: (x) -> (y) : (1 <= x <= 7, x = 2y) These can be thought of as points in the polyhedron: (x, y) : (1 <= x <= 7, x = 2y) This relation contains the pairs (2, 1), (4, 2), and (6, 3).
Since IntegerRelation makes a distinction between dimensions, VarKind::Range and VarKind::Domain should be used to refer to dimension variables.
Definition at line 64 of file IntegerRelation.h.
|
strong |
All derived classes of IntegerRelation.
Enumerator | |
---|---|
IntegerRelation | |
IntegerPolyhedron | |
FlatLinearConstraints | |
FlatLinearValueConstraints | |
FlatAffineValueConstraints | |
FlatAffineRelation |
Definition at line 67 of file IntegerRelation.h.
|
inline |
Constructs a relation reserving memory for the specified number of constraints and variables.
Definition at line 78 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getNumVars(), and space.
Referenced by getConstantBound(), and getUniverse().
|
inlineexplicit |
Constructs a relation with the specified number of dimensions and symbols.
Definition at line 89 of file IntegerRelation.h.
|
virtualdefault |
void IntegerRelation::addBound | ( | BoundType | type, |
ArrayRef< DynamicAPInt > | expr, | ||
const DynamicAPInt & | value | ||
) |
Adds a constant bound for the specified expression.
Definition at line 1498 of file IntegerRelation.cpp.
References mlir::presburger::Matrix< T >::appendExtraRow(), mlir::presburger::EQ, getNumCols(), mlir::presburger::Matrix< T >::getNumRows(), inequalities, and mlir::presburger::LB.
|
inline |
Definition at line 458 of file IntegerRelation.h.
References addBound(), and mlir::presburger::getDynamicAPIntVec().
void IntegerRelation::addBound | ( | BoundType | type, |
unsigned | pos, | ||
const DynamicAPInt & | value | ||
) |
Adds a constant bound for the specified variable.
Definition at line 1483 of file IntegerRelation.cpp.
References mlir::presburger::Matrix< T >::appendExtraRow(), mlir::presburger::EQ, equalities, getNumCols(), inequalities, and mlir::presburger::LB.
Referenced by addBound(), and mlir::presburger::MultiAffineFunction::getLexSet().
|
inline |
Definition at line 451 of file IntegerRelation.h.
References addBound().
void IntegerRelation::addEquality | ( | ArrayRef< DynamicAPInt > | eq | ) |
Adds an equality from the coefficients specified in eq
.
Definition at line 368 of file IntegerRelation.cpp.
References mlir::presburger::Matrix< T >::appendExtraRow(), equalities, and getNumCols().
Referenced by addEquality(), addOrderingConstraints(), mlir::presburger::LinearTransform::applyTo(), computeDirectionVector(), mlir::FlatLinearValueConstraints::FlatLinearValueConstraints(), mlir::presburger::MultiAffineFunction::getAsRelation(), getEmpty(), mlir::presburger::MultiAffineFunction::getLexSet(), mlir::affine::getRelationFromMap(), and removeDuplicateConstraints().
|
inline |
Definition at line 318 of file IntegerRelation.h.
References addEquality(), and mlir::presburger::getDynamicAPIntVec().
void IntegerRelation::addInequality | ( | ArrayRef< DynamicAPInt > | inEq | ) |
Adds an inequality (>= 0) from the coefficients specified in inEq
.
Definition at line 375 of file IntegerRelation.cpp.
References mlir::presburger::Matrix< T >::appendExtraRow(), getNumCols(), and inequalities.
Referenced by addDivisionConstraints(), addInequality(), addLocalFloorDiv(), mlir::scf::addLoopRangeConstraints(), addOrderingConstraints(), mlir::presburger::LinearTransform::applyTo(), mlir::ValueBoundsConstraintSet::comparePos(), mlir::presburger::detail::computePolytopeGeneratingFunction(), mlir::presburger::SymbolicLexSimplex::computeSymbolicIntegerLexMin(), mlir::FlatLinearValueConstraints::FlatLinearValueConstraints(), mlir::presburger::detail::getDual(), mlir::presburger::MultiAffineFunction::getLexSet(), getSetDifference(), mlir::presburger::IntegerPolyhedron::IntegerPolyhedron(), mlir::scf::rewritePeeledMinMaxOp(), mlir::affine::simplifyConstrainedMinMaxOp(), and unionBoundingBox().
|
inline |
Definition at line 313 of file IntegerRelation.h.
References addInequality(), and mlir::presburger::getDynamicAPIntVec().
void IntegerRelation::addLocalFloorDiv | ( | ArrayRef< DynamicAPInt > | dividend, |
const DynamicAPInt & | divisor | ||
) |
Adds a new local variable as the floordiv of an affine function of other variables, the coefficients of which are provided in dividend
and with respect to a positive constant divisor
.
Adds a new local variable as the floordiv of an affine function of other variables, the coefficients of which are provided in 'dividend' and with respect to a positive constant 'divisor'.
Two constraints are added to the system to capture equivalence with the floordiv: q = dividend floordiv c <=> c*q <= dividend <= c*q + c - 1.
Two constraints are added to the system to capture equivalence with the floordiv. q = expr floordiv c <=> c*q <= expr <= c*q + c - 1.
Definition at line 1514 of file IntegerRelation.cpp.
References addInequality(), appendVar(), mlir::presburger::getDivLowerBound(), mlir::presburger::getDivUpperBound(), getNumCols(), and mlir::presburger::Local.
Referenced by addLocalFloorDiv().
|
inline |
Definition at line 469 of file IntegerRelation.h.
References addLocalFloorDiv(), and mlir::presburger::getDynamicAPIntVec().
void IntegerRelation::append | ( | const IntegerRelation & | other | ) |
Appends constraints from other
into this
.
This is equivalent to an intersection with no simplification of any sort attempted.
Definition at line 85 of file IntegerRelation.cpp.
References equalities, getNumEqualities(), getNumInequalities(), mlir::presburger::Matrix< T >::getNumRows(), getSpace(), inequalities, mlir::presburger::PresburgerSpace::isEqual(), mlir::presburger::Matrix< T >::reserveRows(), and space.
Referenced by mlir::affine::FlatAffineRelation::compose(), mlir::FlatLinearValueConstraints::FlatLinearValueConstraints(), mlir::affine::MemRefAccess::getAccessRelation(), intersect(), intersectDomain(), intersectRange(), mergeAndCompose(), and unionBoundingBox().
unsigned IntegerRelation::appendVar | ( | VarKind | kind, |
unsigned | num = 1 |
||
) |
Append num
variables of the specified kind after the last variable of that kind.
The coefficient columns corresponding to the added variables are initialized to zero. Return the absolute column position (i.e., not relative to the kind of variable) of the first appended variable.
Definition at line 363 of file IntegerRelation.cpp.
References getNumVarKind(), and insertVar().
Referenced by addLocalFloorDiv(), mlir::FlatLinearConstraints::appendDimVar(), mlir::FlatLinearConstraints::appendLocalVar(), mlir::FlatLinearConstraints::appendSymbolVar(), compose(), mlir::affine::MemRefAccess::getAccessRelation(), mlir::ValueBoundsConstraintSet::insert(), intersectDomain(), and intersectRange().
void IntegerRelation::applyDomain | ( | const IntegerRelation & | rel | ) |
Given a relation rel
, apply the relation to the domain of this relation.
R1: i -> j : (0 <= i < 2, j = i) R2: i -> k : (k = i floordiv 2) R3: k -> j : (0 <= k < 1, 2k <= j <= 2k + 1)
R1 = {(0, 0), (1, 1)}. R2 maps both 0 and 1 to 0. So R3 = {(0, 0), (0, 1)}.
Formally, R1.applyDomain(R2) = R2.inverse().compose(R1).
Definition at line 2495 of file IntegerRelation.cpp.
void IntegerRelation::applyRange | ( | const IntegerRelation & | rel | ) |
Given a relation rel
, apply the relation to the range of this relation.
Formally, R1.applyRange(R2) is the same as R1.compose(R2) but we provide this for uniformity with applyDomain
.
Definition at line 2501 of file IntegerRelation.cpp.
References compose().
|
inline |
Definition at line 185 of file IntegerRelation.h.
References equalities.
|
inline |
Returns the value at the specified equality row and column.
Definition at line 177 of file IntegerRelation.h.
References equalities.
Referenced by computeConstantLowerOrUpperBound(), constantFoldVar(), eliminateFromConstraint(), findConstraintWithNonZeroAt(), findSymbolicIntegerLexMax(), gaussianEliminate(), getBoundedDirections(), getConstantBoundOnDimSize(), getDivRepr(), hasInvalidConstraint(), isEmptyByGCDTest(), isObviouslyEqual(), print(), and removeRedundantLocalVars().
|
inline |
The same, but casts to int64_t.
This is unsafe and will assert-fail if the value does not fit in an int64_t.
Definition at line 182 of file IntegerRelation.h.
References equalities.
|
inline |
Definition at line 196 of file IntegerRelation.h.
References inequalities.
|
inline |
Returns the value at the specified inequality row and column.
Definition at line 188 of file IntegerRelation.h.
References inequalities.
Referenced by createSeparationCondition(), eliminateFromConstraint(), findConstraintWithNonZeroAt(), findSymbolicIntegerLexMax(), gcdTightenInequalities(), getBoundedDirections(), getDivRepr(), mlir::presburger::detail::getDual(), hasInvalidConstraint(), isObviouslyEqual(), print(), and removeDuplicateConstraints().
|
inline |
The same, but casts to int64_t.
This is unsafe and will assert-fail if the value does not fit in an int64_t.
Definition at line 193 of file IntegerRelation.h.
References inequalities.
|
inlinestatic |
Definition at line 114 of file IntegerRelation.h.
|
virtual |
Replaces the contents of this IntegerRelation with other
.
Definition at line 563 of file IntegerRelation.cpp.
Referenced by getFlattenedAffineExprs().
void IntegerRelation::clearConstraints | ( | ) |
Removes all equalities and inequalities.
Definition at line 479 of file IntegerRelation.cpp.
References equalities, inequalities, and mlir::presburger::Matrix< T >::resizeVertically().
Referenced by unionBoundingBox().
std::unique_ptr< IntegerRelation > IntegerRelation::clone | ( | ) | const |
Definition at line 50 of file IntegerRelation.cpp.
void IntegerRelation::compose | ( | const IntegerRelation & | rel | ) |
Let the relation this
be R1, and the relation rel
be R2.
Modifies R1 to be the composition of R1 and R2: R1;R2.
Formally, if R1: A -> B, and R2: B -> C, then this function returns a relation R3: A -> C such that a point (a, c) belongs to R3 iff there exists b such that (a, b) is in R1 and, (b, c) is in R2.
Definition at line 2469 of file IntegerRelation.cpp.
References appendVar(), convertVarKind(), mlir::presburger::Domain, getDomainSet(), getNumRangeVars(), getRangeSet(), getSpace(), intersectRange(), mlir::presburger::Local, and mlir::presburger::Range.
Referenced by applyDomain(), applyRange(), and mlir::presburger::PresburgerRelation::compose().
|
protected |
Returns the constant lower bound if isLower is true, and the upper bound if isLower is false.
Definition at line 1721 of file IntegerRelation.cpp.
References atEq(), findEqualityToConstant(), getNumCols(), getNumVars(), and projectOut().
|
inlineprotected |
The same, but casts to int64_t.
This is unsafe and will assert-fail if the value does not fit in an int64_t.
Definition at line 773 of file IntegerRelation.h.
PresburgerRelation IntegerRelation::computeReprWithOnlyDivLocals | ( | ) | const |
Compute an equivalent representation of the same set, such that all local vars in all disjuncts have division representations.
This representation may involve local vars that correspond to divisions, and may also be a union of convex disjuncts.
Definition at line 223 of file IntegerRelation.cpp.
References mlir::presburger::SymbolicLexSimplex::computeSymbolicIntegerLexMin(), copy(), mlir::presburger::PWMAFunction::getDomain(), getNumLocalVars(), mlir::presburger::PresburgerSpace::getSetSpace(), getSpace(), mlir::presburger::SymbolicLexOpt::lexopt, mlir::presburger::Local, mlir::presburger::PresburgerSpace::removeVarRange(), mlir::presburger::PresburgerRelation::setSpace(), space, mlir::presburger::SymbolicLexOpt::unboundedDomain, and mlir::presburger::PresburgerSet::unionSet().
std::optional< DynamicAPInt > IntegerRelation::computeVolume | ( | ) | const |
Compute an overapproximation of the number of integer points in the relation.
Symbol vars currently not supported. If the computed overapproximation is infinite, an empty optional is returned.
Definition at line 1241 of file IntegerRelation.cpp.
References mlir::presburger::Simplex::computeIntegerBounds(), getNumDimAndSymbolVars(), getNumSymbolVars(), getNumVars(), mlir::presburger::SimplexBase::isEmpty(), max(), and min().
LogicalResult IntegerRelation::constantFoldVar | ( | unsigned | pos | ) |
Tries to fold the specified variable to a constant using a trivial equality detection; if successful, the constant is substituted for the variable everywhere in the constraint system and then removed from the system.
Definition at line 1559 of file IntegerRelation.cpp.
References atEq(), findEqualityToConstant(), getNumCols(), getNumVars(), and setAndEliminate().
Referenced by constantFoldVarRange().
void IntegerRelation::constantFoldVarRange | ( | unsigned | pos, |
unsigned | num | ||
) |
This method calls constantFoldVar
for the specified range of variables, num
variables starting at position pos
.
Definition at line 1572 of file IntegerRelation.cpp.
References constantFoldVar().
bool IntegerRelation::containsPoint | ( | ArrayRef< DynamicAPInt > | point | ) | const |
Returns true if the given point satisfies the constraints, or false otherwise.
A point satisfies an equality iff the value of the equality at the expression is zero, and it satisfies an inequality iff the value of the inequality at that point is non-negative.
Takes the values of all vars including locals.
Definition at line 991 of file IntegerRelation.cpp.
References getEquality(), getInequality(), getNumEqualities(), getNumInequalities(), and valueAt().
Referenced by containsPoint().
|
inline |
Definition at line 426 of file IntegerRelation.h.
References containsPoint(), and mlir::presburger::getDynamicAPIntVec().
std::optional< SmallVector< DynamicAPInt, 8 > > IntegerRelation::containsPointNoLocal | ( | ArrayRef< DynamicAPInt > | point | ) | const |
Given the values of non-local vars, return a satisfying assignment to the local if one exists, or an empty optional otherwise.
Just substitute the values given and check if an integer sample exists for the local vars.
TODO: this could be made more efficient by handling divisions separately. Instead of finding an integer sample over all the locals, we can first compute the values of the locals that have division representations and only use the integer emptiness check for the locals that don't have this. Handling this correctly requires ordering the divs, though.
Definition at line 1012 of file IntegerRelation.cpp.
References copy(), getNumLocalVars(), getNumVars(), getVarKindOffset(), and mlir::presburger::Local.
Referenced by mlir::presburger::PresburgerRelation::containsPoint(), and containsPointNoLocal().
|
inline |
Definition at line 434 of file IntegerRelation.h.
References containsPointNoLocal(), and mlir::presburger::getDynamicAPIntVec().
|
inline |
Definition at line 609 of file IntegerRelation.h.
References convertVarKind(), and mlir::presburger::Local.
Referenced by mlir::affine::FlatAffineRelation::compose(), mlir::affine::FlatAffineRelation::getDomainSet(), mlir::affine::FlatAffineRelation::getRangeSet(), and mergeAndCompose().
|
inline |
Definition at line 604 of file IntegerRelation.h.
References convertVarKind(), and getNumVarKind().
void IntegerRelation::convertVarKind | ( | VarKind | srcKind, |
unsigned | varStart, | ||
unsigned | varLimit, | ||
VarKind | dstKind, | ||
unsigned | pos | ||
) |
Converts variables of kind srcKind in the range [varStart, varLimit) to variables of kind dstKind.
If pos
is given, the variables are placed at position pos
of dstKind, otherwise they are placed after all the other variables of kind dstKind. The internal ordering among the moved variables is preserved.
Definition at line 1462 of file IntegerRelation.cpp.
References mlir::presburger::PresburgerSpace::convertVarKind(), equalities, getNumVarKind(), getVarKindOffset(), inequalities, mlir::presburger::Matrix< T >::moveColumns(), and space.
Referenced by mlir::affine::checkMemrefAccessDependence(), compose(), convertToLocal(), convertVarKind(), mlir::affine::MemRefAccess::getAccessRelation(), getDomainSet(), getRangeSet(), and inverse().
void IntegerRelation::dump | ( | ) | const |
Definition at line 2619 of file IntegerRelation.cpp.
References print().
Referenced by mlir::affine::checkMemrefAccessDependence(), mlir::ValueBoundsConstraintSet::dump(), fourierMotzkinEliminate(), and mlir::affine::ComputationSliceState::isSliceValid().
|
virtual |
Eliminate the posB^th
local variable, replacing every instance of it with the posA^th
local variable.
This should be used when the two local variables are known to always take the same values.
Definition at line 1297 of file IntegerRelation.cpp.
References mlir::presburger::Matrix< T >::addToColumn(), equalities, getNumLocalVars(), getVarKindOffset(), inequalities, mlir::presburger::Local, and removeVar().
Referenced by mergeLocalVars(), and removeDuplicateDivs().
|
protected |
Searches for a constraint with a non-zero coefficient at colIdx
in equality (isEq=true) or inequality (isEq=false) constraints.
Returns true and sets row found in search in rowIdx
, false otherwise.
Definition at line 570 of file IntegerRelation.cpp.
References atEq(), atIneq(), getNumCols(), getNumEqualities(), and getNumInequalities().
Referenced by gaussianEliminate(), gaussianEliminateVars(), and isColZero().
MaybeOptimum< SmallVector< DynamicAPInt, 8 > > IntegerRelation::findIntegerLexMin | ( | ) | const |
Same as above, but returns lexicographically minimal integer point.
Note: this should be used only when the lexmin is really required. For a generic integer sampling operation, findIntegerSample is more robust and should be preferred. Note that Domain is minimized first, then range.
Definition at line 162 of file IntegerRelation.cpp.
References mlir::presburger::LexSimplex::findIntegerLexMin(), mlir::presburger::MaybeOptimum< T >::getKind(), getNumDimAndSymbolVars(), getNumSymbolVars(), getNumVars(), and mlir::presburger::MaybeOptimum< T >::isBounded().
std::optional< SmallVector< DynamicAPInt, 8 > > IntegerRelation::findIntegerSample | ( | ) | const |
Find an integer sample point satisfying the constraints using a branch and bound algorithm with generalized basis reduction, with some additional processing using Simplex for unbounded sets.
Let this set be S.
Returns an integer sample point if one exists, or an empty Optional otherwise. The returned value also includes values of local ids.
If S is bounded then we directly call into the GBR sampling algorithm. Otherwise, there are some unbounded directions, i.e., vectors v such that S extends to infinity along v or -v. In this case we use an algorithm described in the integer set library (isl) manual and used by the isl_set_sample function in that library. The algorithm is:
1) Apply a unimodular transform T to S to obtain S*T, such that all dimensions in which S*T is bounded lie in the linear span of a prefix of the dimensions.
2) Construct a set B by removing all constraints that involve the unbounded dimensions and then deleting the unbounded dimensions. Note that B is a Bounded set.
3) Try to obtain a sample from B using the GBR sampling algorithm. If no sample is found, return that S is empty.
4) Otherwise, substitute the obtained sample into S*T to obtain a set C. C is a full-dimensional Cone and always contains a sample.
5) Obtain an integer sample from C.
6) Return T*v, where v is the concatenation of the samples from B and C.
The following is a sketch of a proof that a) If the algorithm returns empty, then S is empty. b) If the algorithm returns a sample, it is a valid sample in S.
The algorithm returns empty only if B is empty, in which case S*T is certainly empty since B was obtained by removing constraints and then deleting unconstrained dimensions from S*T. Since T is unimodular, a vector v is in S*T iff T*v is in S. So in this case, since S*T is empty, S is empty too.
Otherwise, the algorithm substitutes the sample from B into S*T. All the constraints of S*T that did not involve unbounded dimensions are satisfied by this substitution. All dimensions in the linear span of the dimensions outside the prefix are unbounded in S*T (step 1). Substituting values for the bounded dimensions cannot make these dimensions bounded, and these are the only remaining dimensions in C, so C is unbounded along every vector (in the positive or negative direction, or both). C is hence a full-dimensional cone and therefore always contains an integer point.
Concatenating the samples from B and C gives a sample v in S*T, so the returned sample T*v is a sample in S.
Definition at line 871 of file IntegerRelation.cpp.
References mlir::presburger::Simplex::findIntegerSample(), mlir::presburger::SimplexBase::isEmpty(), isEmptyByGCDTest(), and mlir::presburger::Simplex::isUnbounded().
Referenced by isIntegerEmpty().
MaybeOptimum< SmallVector< Fraction, 8 > > IntegerRelation::findRationalLexMin | ( | ) | const |
Get the lexicographically minimum rational point satisfying the constraints.
Returns an empty optional if the relation is empty or if the lexmin is unbounded. Symbols are not supported and will result in assert-failure. Note that Domain is minimized first, then range.
Definition at line 142 of file IntegerRelation.cpp.
References mlir::presburger::LexSimplex::findRationalLexMin(), getNumDimAndSymbolVars(), getNumSymbolVars(), getNumVars(), and mlir::presburger::MaybeOptimum< T >::isBounded().
SymbolicLexOpt IntegerRelation::findSymbolicIntegerLexMax | ( | ) | const |
Same as findSymbolicIntegerLexMin but produces lexmax instead of lexmin.
findSymbolicIntegerLexMax is implemented using findSymbolicIntegerLexMin as follows:
this
relation with the sign of each dimension variable in range flipped;this
relation;Definition at line 318 of file IntegerRelation.cpp.
References atEq(), atIneq(), mlir::presburger::PWMAFunction::Piece::domain, findSymbolicIntegerLexMin(), mlir::presburger::PWMAFunction::getAllPieces(), getNumDomainVars(), getNumEqualities(), getNumInequalities(), getNumRangeVars(), mlir::presburger::Matrix< T >::getNumRows(), mlir::presburger::PWMAFunction::getSpace(), mlir::presburger::SymbolicLexOpt::lexopt, mlir::presburger::Matrix< T >::negateRow(), and mlir::presburger::SymbolicLexOpt::unboundedDomain.
SymbolicLexOpt IntegerRelation::findSymbolicIntegerLexMin | ( | ) | const |
Compute the symbolic integer lexmin of the relation.
This finds, for every assignment to the symbols and domain, the lexicographically minimum value attained by the range.
For example, the symbolic lexmin of the set
(x, y)[a, b, c] : (a <= x, b <= x, x <= c)
can be written as
x = a if b <= a, a <= c x = b if a < b, b <= c
This function is stored in the lexopt
function in the result. Some assignments to the symbols might make the set empty. Such points are not part of the function's domain. In the above example, this happens when max(a, b) > c.
For some values of the symbols, the lexmin may be unbounded. SymbolicLexOpt
stores these parts of the symbolic domain in a separate PresburgerSet
, unboundedDomain
.
Definition at line 284 of file IntegerRelation.cpp.
References mlir::presburger::SymbolicLexSimplex::computeSymbolicIntegerLexMin(), mlir::presburger::Domain, getNumDomainVars(), getNumLocalVars(), mlir::presburger::PWMAFunction::getNumOutputs(), getNumSymbolVars(), getNumVars(), mlir::presburger::PresburgerSpace::getSetSpace(), getVarKindEnd(), getVarKindOffset(), mlir::presburger::SymbolicLexOpt::lexopt, mlir::presburger::PWMAFunction::removeOutputs(), and mlir::presburger::Symbol.
Referenced by findSymbolicIntegerLexMax().
|
protectedvirtual |
Eliminates the variable at the specified position using Fourier-Motzkin variable elimination, but uses Gaussian elimination if there is an equality involving that variable.
If the result of the elimination is integer exact, *isResultIntegerExact
is set to true. If darkShadow
is set to true, a potential under approximation (subset) of the rational shadow / exact integer shadow is computed.
Create the new system which has one variable less.
Definition at line 1940 of file IntegerRelation.cpp.
References dump(), getNumVars(), and hasConsistentState().
Referenced by isEmpty(), projectOut(), and mlir::FlatLinearValueConstraints::projectOut().
|
protected |
Perform a Gaussian elimination operation to reduce all equations to standard form.
Returns whether the constraint system was modified.
Definition at line 1125 of file IntegerRelation.cpp.
References atEq(), eliminateFromConstraint(), equalities, findConstraintWithNonZeroAt(), gcdTightenInequalities(), getEmpty(), getNumEqualities(), getNumInequalities(), getNumVars(), getSpace(), inequalities, mlir::presburger::IntMatrix::normalizeRow(), removeEqualityRange(), and mlir::presburger::Matrix< T >::swapRows().
Referenced by simplify().
|
inlineprotected |
Eliminates a single variable at position
from equality and inequality constraints.
Returns success
if the variable was eliminated, and failure
otherwise.
Definition at line 781 of file IntegerRelation.h.
References gaussianEliminateVars().
|
protected |
Eliminates variables from equality and inequality constraints in column range [posStart, posLimit).
Returns the number of variables eliminated.
Definition at line 1077 of file IntegerRelation.cpp.
References eliminateFromConstraint(), equalities, findConstraintWithNonZeroAt(), gcdTightenInequalities(), getNumEqualities(), getNumInequalities(), getNumVars(), hasConsistentState(), inequalities, mlir::presburger::IntMatrix::normalizeRow(), removeEquality(), and removeVarRange().
Referenced by gaussianEliminateVar(), isEmpty(), and projectOut().
|
protected |
Tightens inequalities given that we are dealing with integer spaces.
This is similar to the GCD test but applied to inequalities. The constant term can be reduced to the preceding multiple of the GCD of the coefficients, i.e., 64*i - 100 >= 0 => 64*i - 128 >= 0 (since 'i' is an integer). This is a fast method (linear in the number of coefficients).
This is analogous to the GCD test but applied to inequalities. The constant term can be reduced to the preceding multiple of the GCD of the coefficients, i.e., 64*i - 100 >= 0 => 64*i - 128 >= 0 (since 'i' is an integer). This is a fast method - linear in the number of coefficients.
Definition at line 1065 of file IntegerRelation.cpp.
References atIneq(), getNumCols(), getNumInequalities(), inequalities, and mlir::presburger::IntMatrix::normalizeRow().
Referenced by gaussianEliminate(), gaussianEliminateVars(), projectOut(), removeRedundantConstraints(), and removeTrivialRedundancy().
IntMatrix IntegerRelation::getBoundedDirections | ( | ) | const |
Returns a matrix where each row is a vector along which the polytope is bounded.
The span of the returned vectors is guaranteed to contain all such vectors. The returned vectors are NOT guaranteed to be linearly independent. This function should not be called on empty sets.
Definition at line 782 of file IntegerRelation.cpp.
References atEq(), atIneq(), getNumCols(), getNumEqualities(), getNumInequalities(), mlir::presburger::Simplex::isBoundedAlongConstraint(), and mlir::presburger::SimplexBase::isEmpty().
std::optional< DynamicAPInt > IntegerRelation::getConstantBound | ( | BoundType | type, |
unsigned | pos | ||
) | const |
Returns the constant bound for the pos^th variable if there is one; std::nullopt otherwise.
Definition at line 1778 of file IntegerRelation.cpp.
References mlir::presburger::EQ, IntegerRelation(), mlir::presburger::LB, and mlir::presburger::UB.
Referenced by getConstantBound64(), and unionBoundingBox().
|
inline |
The same, but casts to int64_t.
This is unsafe and will assert-fail if the value does not fit in an int64_t.
Definition at line 558 of file IntegerRelation.h.
References getConstantBound().
Referenced by mlir::ValueBoundsConstraintSet::computeConstantBound(), computeDirectionVector(), mlir::affine::normalizeMemRefType(), and mlir::affine::simplifyConstrainedMinMaxOp().
std::optional< DynamicAPInt > IntegerRelation::getConstantBoundOnDimSize | ( | unsigned | pos, |
SmallVectorImpl< DynamicAPInt > * | lb = nullptr , |
||
DynamicAPInt * | boundFloorDivisor = nullptr , |
||
SmallVectorImpl< DynamicAPInt > * | ub = nullptr , |
||
unsigned * | minLbPos = nullptr , |
||
unsigned * | minUbPos = nullptr |
||
) | const |
Returns the smallest known constant bound for the extent of the specified variable (pos^th), i.e., the smallest known constant that is greater than or equal to 'exclusive upper bound' - 'lower bound' of the variable.
Returns a non-negative constant bound on the extent (upper bound - lower bound) of the specified variable if it is found to be a constant; returns std::nullopt if it's not a constant.
This constant bound is guaranteed to be non-negative. Returns std::nullopt if it's not a constant. This method employs trivial (low complexity / cost) checks and detection. Symbolic variables are treated specially, i.e., it looks for constant differences between affine expressions involving only the symbolic variables. lb
and ub
(along with the boundFloorDivisor
) are set to represent the lower and upper bound associated with the constant difference: lb
, ub
have the coefficients, and boundFloorDivisor
, their divisor. minLbPos
and minUbPos
if non-null are set to the position of the constant lower bound and upper bound respectively (to the same if they are from an equality). Ex: if the lower bound is [(s0 + s2 - 1) floordiv 32] for a system with three symbolic variables, *lb = [1, 0, 1], lbDivisor = 32. See comments at function definition for examples.
This methods treats symbolic variables specially, i.e., it looks for constant differences between affine expressions involving only the symbolic variables. See comments at function definition for example. 'lb', if provided, is set to the lower bound associated with the constant difference. Note that 'lb' is purely symbolic and thus will contain the coefficients of the symbolic variables and the constant coefficient.
Definition at line 1593 of file IntegerRelation.cpp.
References atEq(), findEqualityToConstant(), getEquality(), getNumDimAndSymbolVars(), getNumDimVars(), and getNumSymbolVars().
Referenced by createFullTiles(), createSeparationCondition(), getConstantBoundOnDimSize64(), and unionBoundingBox().
|
inline |
The same, but casts to int64_t.
This is unsafe and will assert-fail if the value does not fit in an int64_t.
Definition at line 533 of file IntegerRelation.h.
References getConstantBoundOnDimSize(), and mlir::presburger::getInt64Vec().
Referenced by mlir::affine::MemRefRegion::getConstantBoundOnDimSize().
IntegerRelation::CountsSnapshot IntegerRelation::getCounts | ( | ) | const |
Definition at line 199 of file IntegerRelation.cpp.
References getNumEqualities(), getNumInequalities(), and getSpace().
Referenced by mlir::presburger::SymbolicLexSimplex::computeSymbolicIntegerLexMin(), and getSetDifference().
IntegerPolyhedron IntegerRelation::getDomainSet | ( | ) | const |
Return a set corresponding to all points in the domain of the relation.
Definition at line 2346 of file IntegerRelation.cpp.
References convertVarKind(), mlir::presburger::Domain, getNumVarKind(), mlir::presburger::Local, mlir::presburger::Range, and mlir::presburger::SetDim.
Referenced by mlir::affine::checkMemrefAccessDependence(), compose(), and intersectDomain().
|
inlinestatic |
Return an empty system containing an invalid equation 0 = 1.
Definition at line 103 of file IntegerRelation.h.
References addEquality(), mlir::presburger::PresburgerSpace::getNumVars(), and space.
Referenced by gaussianEliminate(), and removeDuplicateConstraints().
|
inline |
Definition at line 232 of file IntegerRelation.h.
References equalities, and mlir::presburger::Matrix< T >::getRow().
Referenced by mlir::presburger::LinearTransform::applyTo(), containsPoint(), getConstantBoundOnDimSize(), getIneqCoeffsFromIdx(), mlir::presburger::SimplexBase::intersectIntegerRelation(), mlir::presburger::Simplex::isRationalSubsetOf(), removeConstraintsInvolvingVarRange(), and removeTrivialEqualities().
|
inline |
The same, but casts to int64_t.
This is unsafe and will assert-fail if the value does not fit in an int64_t.
Definition at line 240 of file IntegerRelation.h.
References equalities, mlir::presburger::getInt64Vec(), and mlir::presburger::Matrix< T >::getRow().
ArrayRef< Identifier > IntegerRelation::getIds | ( | VarKind | kind | ) |
Get the identifiers for the variables of specified varKind.
Calls resetIds on the relations space if identifiers are not enabled.
Definition at line 79 of file IntegerRelation.cpp.
References mlir::presburger::PresburgerSpace::getIds(), mlir::presburger::PresburgerSpace::isUsingIds(), mlir::presburger::PresburgerSpace::resetIds(), and space.
Referenced by mlir::affine::MemRefAccess::getAccessRelation().
|
inline |
Definition at line 247 of file IntegerRelation.h.
References inequalities.
Referenced by mlir::presburger::detail::computePolytopeGeneratingFunction(), and mlir::presburger::detail::computeUnimodularConeGeneratingFunction().
|
inline |
Definition at line 235 of file IntegerRelation.h.
References mlir::presburger::Matrix< T >::getRow(), and inequalities.
Referenced by mlir::presburger::LinearTransform::applyTo(), containsPoint(), createSeparationCondition(), getIneqCoeffsFromIdx(), getSetDifference(), mlir::presburger::SimplexBase::intersectIntegerRelation(), isFullDim(), mlir::presburger::Simplex::isRationalSubsetOf(), removeConstraintsInvolvingVarRange(), and removeDuplicateConstraints().
|
inline |
Definition at line 243 of file IntegerRelation.h.
References mlir::presburger::getInt64Vec(), mlir::presburger::Matrix< T >::getRow(), and inequalities.
|
inlinevirtual |
Return the kind of this IntegerRelation.
Reimplemented in mlir::affine::FlatAffineRelation, mlir::affine::FlatAffineValueConstraints, mlir::presburger::IntegerPolyhedron, mlir::FlatLinearValueConstraints, and mlir::FlatLinearConstraints.
Definition at line 112 of file IntegerRelation.h.
References IntegerRelation.
Referenced by mlir::presburger::IntegerPolyhedron::classof().
DivisionRepr IntegerRelation::getLocalReprs | ( | std::vector< MaybeLocalRepr > * | repr = nullptr | ) | const |
Returns a DivisonRepr
representing the division representation of local variables in the constraint system.
If repr
is not nullptr
, the equality and pairs of inequality constraints identified by their position indices using which an explicit representation for each local variable can be computed are set in repr
in the form of a MaybeLocalRepr
struct. If no such inequality pair/equality can be found, the kind attribute in MaybeLocalRepr
is set to None.
Definition at line 1023 of file IntegerRelation.cpp.
References mlir::presburger::DivisionRepr::clearRepr(), mlir::presburger::computeSingleVarRepr(), mlir::presburger::DivisionRepr::getDenom(), mlir::presburger::DivisionRepr::getDividend(), getNumDimAndSymbolVars(), getNumLocalVars(), getNumVars(), getVarKindOffset(), and mlir::presburger::Local.
Referenced by mlir::getMultiAffineFunctionFromMap(), getSetDifference(), hasOnlyDivLocals(), mlir::presburger::mergeLocalVars(), and removeDuplicateDivs().
void IntegerRelation::getLowerAndUpperBoundIndices | ( | unsigned | pos, |
SmallVectorImpl< unsigned > * | lbIndices, | ||
SmallVectorImpl< unsigned > * | ubIndices, | ||
SmallVectorImpl< unsigned > * | eqIndices = nullptr , |
||
unsigned | offset = 0 , |
||
unsigned | num = 0 |
||
) | const |
Gather positions of all lower and upper bounds of the variable at pos
, and optionally any equalities on it.
Gather all lower and upper bounds of the variable at pos
, and optionally any equalities on it.
In addition, the bounds are to be independent of variables in position range [offset
, offset
+ num
).
Definition at line 487 of file IntegerRelation.cpp.
References getNumCols(), and getNumVars().
Referenced by mlir::presburger::computeSingleVarRepr(), and createSeparationCondition().
|
inline |
Returns the number of columns in the constraint system.
Definition at line 216 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getNumVars(), and space.
Referenced by addBound(), addEquality(), addInequality(), addLocalFloorDiv(), mlir::scf::addLoopRangeConstraints(), addOrderingConstraints(), mlir::ValueBoundsConstraintSet::comparePos(), computeConstantLowerOrUpperBound(), computeDirectionVector(), mlir::presburger::computeSingleVarRepr(), constantFoldVar(), createPrivateMemRef(), createSeparationCondition(), detectAsFloorDiv(), eliminateFromConstraint(), findConstraintWithNonZeroAt(), mlir::FlatLinearValueConstraints::FlatLinearValueConstraints(), gcdTightenInequalities(), generateCopy(), mlir::presburger::MultiAffineFunction::getAsRelation(), getBoundedDirections(), getDivRepr(), mlir::presburger::detail::getDual(), getLowerAndUpperBoundIndices(), mlir::affine::getRelationFromMap(), hasInvalidConstraint(), isEmptyByGCDTest(), isHyperRectangular(), isObviouslyEqual(), print(), projectOut(), removeDuplicateConstraints(), setAndEliminate(), mlir::affine::simplifyConstrainedMinMaxOp(), and unionBoundingBox().
|
inline |
Definition at line 200 of file IntegerRelation.h.
References getNumEqualities(), and getNumInequalities().
Referenced by isEmpty(), and printSpace().
|
inline |
Definition at line 210 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getNumDimAndSymbolVars(), and space.
Referenced by areVarsUnique(), mlir::ValueBoundsConstraintSet::computeBound(), mlir::affine::computeSliceUnion(), computeVolume(), createFullTiles(), createSeparationCondition(), findIntegerLexMin(), findRationalLexMin(), mlir::FlatLinearValueConstraints::FlatLinearValueConstraints(), mlir::affine::getComputationSliceState(), getConstantBoundOnDimSize(), getLocalReprs(), mlir::FlatLinearValueConstraints::getMaybeValues(), mlir::ValueBoundsConstraintSet::getPosExpr(), mlir::affine::getRelationFromMap(), mlir::FlatLinearValueConstraints::getValue(), mlir::FlatLinearValueConstraints::getValues(), mlir::FlatLinearValueConstraints::hasValue(), mlir::FlatLinearValueConstraints::mergeSymbolVars(), mlir::FlatLinearConstraints::printSpace(), mlir::FlatLinearValueConstraints::printSpace(), removeRedundantLocalVars(), mlir::FlatLinearValueConstraints::setValue(), and mlir::affine::simplifyConstrainedMinMaxOp().
|
inline |
Definition at line 209 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getNumDimVars(), and space.
Referenced by mlir::ValueBoundsConstraintSet::addBound(), mlir::scf::addLoopRangeConstraints(), addMissingLoopIVBounds(), addOrderingConstraints(), mlir::FlatLinearValueConstraints::appendDimVar(), areVarsUnique(), mlir::FlatLinearValueConstraints::computeAlignedMap(), mlir::ValueBoundsConstraintSet::computeBound(), computeDirectionVector(), mlir::affine::computeSliceUnion(), mlir::presburger::SymbolicLexSimplex::computeSymbolicIntegerLexMin(), generateCopy(), mlir::affine::MemRefAccess::getAccessRelation(), getConstantBoundOnDimSize(), getNumCommonLoops(), mlir::ValueBoundsConstraintSet::getPosExpr(), mlir::affine::getRelationFromMap(), mlir::affine::ComputationSliceState::isMaximal(), mergeAndAlignVars(), mlir::FlatLinearValueConstraints::mergeSymbolVars(), turnSymbolIntoDim(), mlir::FlatLinearValueConstraints::unionBoundingBox(), and unionBoundingBox().
|
inline |
Definition at line 204 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getNumDomainVars(), and space.
Referenced by areVarsAligned(), mlir::affine::checkMemrefAccessDependence(), findSymbolicIntegerLexMax(), findSymbolicIntegerLexMin(), intersectRange(), and mergeAndCompose().
|
inline |
Definition at line 218 of file IntegerRelation.h.
References equalities, and mlir::presburger::Matrix< T >::getNumRows().
Referenced by append(), mlir::presburger::LinearTransform::applyTo(), containsPoint(), findConstraintWithNonZeroAt(), findSymbolicIntegerLexMax(), gaussianEliminate(), gaussianEliminateVars(), getBoundedDirections(), getCounts(), getDivRepr(), getIneqCoeffsFromIdx(), getNumConstraints(), getSetDifference(), hasInvalidConstraint(), mlir::presburger::SimplexBase::intersectIntegerRelation(), isEmptyByGCDTest(), isFullDim(), isObviouslyEqual(), mlir::presburger::Simplex::isRationalSubsetOf(), normalizeConstraintsByGCD(), print(), removeConstraintsInvolvingVarRange(), removeRedundantLocalVars(), removeTrivialEqualities(), and truncate().
|
inline |
Definition at line 220 of file IntegerRelation.h.
References mlir::presburger::Matrix< T >::getNumRows(), and inequalities.
Referenced by append(), mlir::presburger::LinearTransform::applyTo(), mlir::ValueBoundsConstraintSet::comparePos(), mlir::presburger::detail::computePolytopeGeneratingFunction(), mlir::presburger::detail::computeUnimodularConeGeneratingFunction(), containsPoint(), findConstraintWithNonZeroAt(), findSymbolicIntegerLexMax(), gaussianEliminate(), gaussianEliminateVars(), gcdTightenInequalities(), getBoundedDirections(), getCounts(), getDivRepr(), mlir::presburger::detail::getDual(), getIneqCoeffsFromIdx(), mlir::presburger::MultiAffineFunction::getLexSet(), getNumConstraints(), getSetDifference(), hasInvalidConstraint(), mlir::presburger::SimplexBase::intersectIntegerRelation(), isFullDim(), isObviouslyEqual(), mlir::presburger::Simplex::isRationalSubsetOf(), normalizeConstraintsByGCD(), print(), removeConstraintsInvolvingVarRange(), removeDuplicateConstraints(), removeRedundantConstraints(), removeRedundantInequalities(), removeRedundantLocalVars(), and truncate().
|
inline |
Definition at line 207 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getNumLocalVars(), and space.
Referenced by addDivisionConstraints(), computeReprWithOnlyDivLocals(), mlir::affine::computeSliceUnion(), containsPointNoLocal(), eliminateRedundantLocalVar(), findSymbolicIntegerLexMin(), mlir::FlatLinearValueConstraints::FlatLinearValueConstraints(), mlir::FlatLinearConstraints::flattenAlignedMapAndMergeLocals(), getLocalReprs(), getSetDifference(), mlir::affine::ComputationSliceState::isSliceValid(), mergeLocalVars(), mlir::presburger::mergeLocalVars(), mlir::affine::normalizeMemRefType(), mlir::FlatLinearValueConstraints::unionBoundingBox(), and unionBoundingBox().
|
inline |
Definition at line 205 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getNumRangeVars(), and space.
Referenced by areVarsAligned(), compose(), mlir::presburger::detail::computePolytopeGeneratingFunction(), findSymbolicIntegerLexMax(), intersectDomain(), and mergeAndCompose().
|
inline |
Definition at line 224 of file IntegerRelation.h.
References equalities, and mlir::presburger::Matrix< T >::getNumReservedRows().
|
inline |
Definition at line 228 of file IntegerRelation.h.
References mlir::presburger::Matrix< T >::getNumReservedRows(), and inequalities.
|
inline |
Definition at line 206 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getNumSymbolVars(), and space.
Referenced by mlir::ValueBoundsConstraintSet::addBound(), mlir::scf::addLoopRangeConstraints(), mlir::FlatLinearValueConstraints::appendSymbolVar(), areVarsAligned(), mlir::FlatLinearValueConstraints::computeAlignedMap(), mlir::presburger::detail::computePolytopeGeneratingFunction(), mlir::presburger::SymbolicLexSimplex::computeSymbolicIntegerLexMin(), computeVolume(), findIntegerLexMin(), findRationalLexMin(), findSymbolicIntegerLexMin(), generateCopy(), getConstantBoundOnDimSize(), mlir::affine::getRelationFromMap(), mlir::affine::ComputationSliceState::isSliceValid(), mlir::FlatLinearValueConstraints::mergeSymbolVars(), mlir::presburger::SymbolicLexSimplex::SymbolicLexSimplex(), turnSymbolIntoDim(), and unionBoundingBox().
|
inline |
Get the number of vars of the specified kind.
Definition at line 250 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getNumVarKind(), and space.
Referenced by appendVar(), mlir::FlatLinearValueConstraints::computeAlignedMap(), convertVarKind(), mlir::FlatLinearValueConstraints::FlatLinearValueConstraints(), getDomainSet(), mlir::FlatLinearValueConstraints::getMaybeValues(), getRangeSet(), insertVar(), inverse(), mergeAndAlignSymbols(), removeVarRange(), and truncateVarKind().
|
inline |
Definition at line 213 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getNumVars(), and space.
Referenced by addDivisionConstraints(), checkExplicitRepresentation(), computeConstantLowerOrUpperBound(), computeDirectionVector(), mlir::presburger::computeSingleVarRepr(), mlir::presburger::detail::computeUnimodularConeGeneratingFunction(), computeVolume(), constantFoldVar(), containsPointNoLocal(), createPrivateMemRef(), detectAsFloorDiv(), findEqualityToConstant(), findIntegerLexMin(), findRationalLexMin(), findSymbolicIntegerLexMin(), mlir::FlatLinearConstraints::FlatLinearConstraints(), fourierMotzkinEliminate(), gaussianEliminate(), gaussianEliminateVars(), generateCopy(), getBestVarToEliminate(), getDivRepr(), getIndependentConstraints(), getLocalReprs(), getLowerAndUpperBoundIndices(), mlir::presburger::SimplexBase::intersectIntegerRelation(), isEmpty(), isFullDim(), mlir::affine::ComputationSliceState::isSliceValid(), mlir::affine::normalizeMemRefType(), projectOut(), removeIndependentConstraints(), removeRedundantLocalVars(), removeVarRange(), setAndEliminate(), setSpaceExceptLocals(), mlir::FlatLinearValueConstraints::setValues(), swapVar(), and mlir::presburger::SymbolicLexSimplex::SymbolicLexSimplex().
IntegerPolyhedron IntegerRelation::getRangeSet | ( | ) | const |
Return a set corresponding to all points in the range of the relation.
Definition at line 2420 of file IntegerRelation.cpp.
References convertVarKind(), mlir::presburger::Domain, getNumVarKind(), and mlir::presburger::Local.
Referenced by compose(), and intersectRange().
|
inline |
Returns a reference to the underlying space.
Definition at line 120 of file IntegerRelation.h.
References space.
Referenced by append(), mlir::presburger::LinearTransform::applyTo(), mlir::affine::FlatAffineRelation::compose(), compose(), computeReprWithOnlyDivLocals(), gaussianEliminate(), mlir::affine::MemRefAccess::getAccessRelation(), getCommonConstraints(), getCounts(), getSetDifference(), intersectDomain(), intersectRange(), isEqual(), isObviouslyEqual(), isSubsetOf(), mergeAndCompose(), mlir::presburger::mergeLocalVars(), removeDuplicateConstraints(), unionBoundingBox(), and mlir::presburger::PresburgerRelation::unionInPlace().
|
inline |
Returns a copy of the space without locals.
Definition at line 145 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getNumDomainVars(), mlir::presburger::PresburgerSpace::getNumRangeVars(), mlir::presburger::PresburgerSpace::getNumSymbolVars(), mlir::presburger::PresburgerSpace::getRelationSpace(), and space.
Referenced by getSetDifference().
|
inlinestatic |
Return a system with no constraints, i.e., one which is satisfied by all points.
Definition at line 98 of file IntegerRelation.h.
References IntegerRelation(), and space.
Referenced by mlir::presburger::PresburgerRelation::complement(), and mlir::presburger::PresburgerRelation::getUniverse().
|
inline |
Return the VarKind of the var at the specified position.
Definition at line 272 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getVarKindAt(), and space.
Referenced by mlir::FlatLinearValueConstraints::getValue(), mlir::FlatLinearValueConstraints::hasValue(), and mlir::FlatLinearValueConstraints::setValue().
|
inline |
Return the index at Which the specified kind of vars ends.
Definition at line 260 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getVarKindEnd(), and space.
Referenced by findSymbolicIntegerLexMin(), mlir::presburger::MultiAffineFunction::getAsRelation(), inverse(), mlir::FlatLinearConstraints::printSpace(), and mlir::FlatLinearValueConstraints::printSpace().
|
inline |
Return the index at which the specified kind of vars starts.
Definition at line 255 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getVarKindOffset(), and space.
Referenced by containsPointNoLocal(), convertVarKind(), eliminateRedundantLocalVar(), findSymbolicIntegerLexMin(), mlir::presburger::MultiAffineFunction::getAsRelation(), getLocalReprs(), getSetDifference(), mlir::FlatLinearValueConstraints::getValue(), mlir::FlatLinearValueConstraints::hasValue(), mergeAndAlignSymbols(), mlir::FlatLinearConstraints::printSpace(), mlir::FlatLinearValueConstraints::printSpace(), removeVarRange(), mlir::FlatLinearValueConstraints::setValue(), and swapVar().
|
inline |
Get the number of elements of the specified kind in the range [varStart, varLimit).
Definition at line 266 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::getVarKindOverlap(), and space.
|
protectedvirtual |
Returns false if the fields corresponding to various variable counts, or equality/inequality buffer sizes aren't consistent; true otherwise.
This is meant to be used within an assert internally.
Definition at line 538 of file IntegerRelation.cpp.
References equalities, mlir::presburger::Matrix< T >::hasConsistentState(), and inequalities.
Referenced by fourierMotzkinEliminate(), gaussianEliminateVars(), hasInvalidConstraint(), isEmptyByGCDTest(), and print().
|
protected |
Checks all rows of equality/inequality constraints for trivial contradictions (for example: 1 == 0, 0 >= 1), which may have surfaced after elimination.
Returns true if an invalid constraint is found; false otherwise.
Definition at line 592 of file IntegerRelation.cpp.
References atEq(), atIneq(), getNumCols(), getNumEqualities(), getNumInequalities(), and hasConsistentState().
Referenced by isEmpty(), and isObviouslyEmpty().
bool IntegerRelation::hasOnlyDivLocals | ( | ) | const |
Check whether all local ids have a division representation.
Definition at line 1387 of file IntegerRelation.cpp.
References getLocalReprs(), and mlir::presburger::DivisionRepr::hasAllReprs().
Referenced by mlir::presburger::PresburgerRelation::hasOnlyDivLocals().
|
virtual |
Insert num
variables of the specified kind at position pos
.
Positions are relative to the kind of variable. The coefficient columns corresponding to the added variables are initialized to zero. Return the absolute column position (i.e., not relative to the kind of variable) of the first added variable.
Reimplemented in mlir::presburger::IntegerPolyhedron, and mlir::FlatLinearValueConstraints.
Definition at line 354 of file IntegerRelation.cpp.
References equalities, getNumVarKind(), inequalities, mlir::presburger::Matrix< T >::insertColumns(), mlir::presburger::PresburgerSpace::insertVar(), and space.
Referenced by appendVar(), mlir::affine::MemRefAccess::getAccessRelation(), mlir::presburger::MultiAffineFunction::getAsRelation(), mlir::presburger::IntegerPolyhedron::insertVar(), mergeAndAlignSymbols(), mergeAndCompose(), and mlir::presburger::mergeLocalVars().
IntegerRelation IntegerRelation::intersect | ( | IntegerRelation | other | ) | const |
Return the intersection of the two relations.
If there are locals, they will be merged.
Definition at line 100 of file IntegerRelation.cpp.
References append(), and mergeLocalVars().
Referenced by mlir::presburger::IntegerPolyhedron::intersect(), and mlir::presburger::PresburgerRelation::intersect().
void IntegerRelation::intersectDomain | ( | const IntegerPolyhedron & | poly | ) |
Intersect the given poly
with the domain in-place.
Formally, let the relation this
be R: A -> B and poly is C, then this operation modifies R to be (A intersection C) -> B.
Definition at line 2433 of file IntegerRelation.cpp.
References append(), appendVar(), getDomainSet(), getNumRangeVars(), getSpace(), inverse(), mergeLocalVars(), and mlir::presburger::Range.
Referenced by mlir::presburger::MultiAffineFunction::isEqual().
void IntegerRelation::intersectRange | ( | const IntegerPolyhedron & | poly | ) |
Intersect the given poly
with the range in-place.
Formally, let the relation this
be R: A -> B and poly is C, then this operation modifies R to be A -> (B intersection C).
Definition at line 2449 of file IntegerRelation.cpp.
References append(), appendVar(), mlir::presburger::Domain, getNumDomainVars(), getRangeSet(), getSpace(), and mergeLocalVars().
Referenced by compose().
void IntegerRelation::inverse | ( | ) |
Invert the relation i.e., swap its domain and range.
Formally, let the relation this
be R: A -> B, then this operation modifies R to be B -> A.
Definition at line 2462 of file IntegerRelation.cpp.
References convertVarKind(), mlir::presburger::Domain, getNumVarKind(), getVarKindEnd(), and mlir::presburger::Range.
Referenced by applyDomain(), mlir::affine::checkMemrefAccessDependence(), and intersectDomain().
|
protected |
Returns true if the pos^th column is all zero for both inequalities and equalities.
Definition at line 2292 of file IntegerRelation.cpp.
References findConstraintWithNonZeroAt().
bool IntegerRelation::isEmpty | ( | ) | const |
Checks for emptiness by performing variable elimination on all variables, running the GCD test on each equality constraint, and checking for invalid constraints.
Returns true if the GCD test fails for any equality, or if any invalid constraints are discovered on any row. Returns false otherwise.
Definition at line 695 of file IntegerRelation.cpp.
References fourierMotzkinEliminate(), gaussianEliminateVars(), getBestVarToEliminate(), getNumConstraints(), getNumVars(), hasInvalidConstraint(), isEmptyByGCDTest(), kExplosionFactor, and removeRedundantLocalVars().
Referenced by mlir::affine::checkMemrefAccessDependence(), mlir::ValueBoundsConstraintSet::comparePos(), mlir::presburger::PresburgerRelation::compose(), mlir::presburger::PresburgerRelation::intersect(), isFullDim(), mlir::affine::simplifyConstrainedMinMaxOp(), and mlir::affine::simplifyIntegerSet().
bool IntegerRelation::isEmptyByGCDTest | ( | ) | const |
Runs the GCD test on all equality constraints.
Returns true if this test fails on any equality. Returns false otherwise. This test can be used to disprove the existence of a solution. If it returns true, no integer solution to the equality constraints can exist.
Definition at line 759 of file IntegerRelation.cpp.
References mlir::presburger::abs(), atEq(), getNumCols(), getNumEqualities(), and hasConsistentState().
Referenced by findIntegerSample(), getSetDifference(), isEmpty(), and isObviouslyEmpty().
bool IntegerRelation::isEqual | ( | const IntegerRelation & | other | ) | const |
Return whether this
and other
are equal.
This is integer-exact and somewhat expensive, since it uses the integer emptiness check (see IntegerRelation::findIntegerSample()).
Definition at line 107 of file IntegerRelation.cpp.
References getSpace(), mlir::presburger::PresburgerSpace::isCompatible(), mlir::presburger::PresburgerRelation::isEqual(), and space.
Referenced by mlir::presburger::MultiAffineFunction::isEqual().
bool IntegerRelation::isFullDim | ( | ) |
Definition at line 2514 of file IntegerRelation.cpp.
References getInequality(), getNumEqualities(), getNumInequalities(), getNumVars(), isEmpty(), mlir::presburger::Simplex::isFlatAlong(), and removeTrivialEqualities().
Referenced by mlir::presburger::PresburgerRelation::isFullDim().
bool IntegerRelation::isHyperRectangular | ( | unsigned | pos, |
unsigned | num | ||
) | const |
Returns true if the set can be trivially detected as being hyper-rectangular on the specified contiguous set of variables.
Definition at line 1798 of file IntegerRelation.cpp.
References getNumCols().
Referenced by checkIfHyperRectangular().
bool IntegerRelation::isIntegerEmpty | ( | ) | const |
Returns true if the set of constraints is found to have no solution, false if a solution exists.
Uses the same algorithm as findIntegerSample
.
Definition at line 823 of file IntegerRelation.cpp.
References findIntegerSample().
Referenced by mlir::presburger::PresburgerRelation::isIntegerEmpty().
bool IntegerRelation::isObviouslyEmpty | ( | ) | const |
Performs GCD checks and invalid constraint checks.
Definition at line 741 of file IntegerRelation.cpp.
References hasInvalidConstraint(), and isEmptyByGCDTest().
Referenced by simplify(), and mlir::presburger::PresburgerRelation::simplify().
bool IntegerRelation::isObviouslyEqual | ( | const IntegerRelation & | other | ) | const |
Perform a quick equality check on this
and other
.
The relations are equal if the check return true, but may or may not be equal if the check returns false. The equality check is performed in a plain manner, by comparing if all the equalities and inequalities in this
and other
are the same.
Definition at line 112 of file IntegerRelation.cpp.
References atEq(), atIneq(), cols, getNumCols(), getNumEqualities(), getNumInequalities(), getSpace(), mlir::presburger::PresburgerSpace::isEqual(), and space.
bool IntegerRelation::isSubsetOf | ( | const IntegerRelation & | other | ) | const |
Return whether this is a subset of the given IntegerRelation.
This is integer-exact and somewhat expensive, since it uses the integer emptiness check (see IntegerRelation::findIntegerSample()).
Definition at line 136 of file IntegerRelation.cpp.
References getSpace(), mlir::presburger::PresburgerSpace::isCompatible(), mlir::presburger::PresburgerRelation::isSubsetOf(), and space.
void IntegerRelation::mergeAndAlignSymbols | ( | IntegerRelation & | other | ) |
Merge and align symbol variables of this
and other
with respect to identifiers.
mergeAndAlignSymbols's implementation can be broken down into two steps:
After this operation the symbol variables of both relations have the same identifiers in the same order.
other
from this. If an identifier from
thisexists in
otherthen we align it. Otherwise, we assume it is a new identifier and insert it into
otherin the same position as
this.
Add identifiers that are in
otherbut not
this to this
. Definition at line 1314 of file IntegerRelation.cpp.
References mlir::presburger::PresburgerSpace::getId(), mlir::presburger::PresburgerSpace::getIds(), getNumVarKind(), getVarKindOffset(), insertVar(), mlir::presburger::PresburgerSpace::isUsingIds(), mlir::presburger::PresburgerSpace::setId(), space, swapVar(), and mlir::presburger::Symbol.
Referenced by mlir::affine::MemRefAccess::getAccessRelation(), and mergeAndCompose().
void IntegerRelation::mergeAndCompose | ( | const IntegerRelation & | other | ) |
Given a relation other: (A -> B)
, this operation merges the symbol and local variables and then takes the composition of other
on this: (B -> C)
.
The resulting relation represents tuples of the form: A -> C
.
Definition at line 2533 of file IntegerRelation.cpp.
References append(), convertToLocal(), mlir::presburger::Domain, mlir::presburger::PresburgerSpace::getId(), getNumDomainVars(), getNumRangeVars(), getSpace(), mlir::presburger::Identifier::hasValue(), insertVar(), mergeAndAlignSymbols(), mergeLocalVars(), mlir::presburger::Range, removeRedundantLocalVars(), removeVarRange(), mlir::presburger::PresburgerSpace::setId(), and space.
Referenced by mlir::affine::checkMemrefAccessDependence().
unsigned IntegerRelation::mergeLocalVars | ( | IntegerRelation & | other | ) |
Adds additional local vars to the sets such that they both have the union of the local vars in each set, without changing the set of points that lie in this
and other
.
Adds additional local ids to the sets such that they both have the union of the local ids in each set, without changing the set of points that lie in this
and other
.
While taking union, if a local var in other
has a division representation which is a duplicate of division representation, of another local var, it is not added to the final union of local vars and is instead merged. The new ordering of local vars is:
[Local vars of this
] [Non-merged local vars of other
]
The relative ordering of local vars is same as before.
After merging, if the i^th
local variable in one set has a known division representation, then the i^th
local variable in the other set either has the same division representation or no known division representation.
The spaces of both relations should be compatible.
Returns the number of non-merged local vars of other
, i.e. the number of locals that have been added to this
.
To detect local ids that always take the same value, each local id is represented as a floordiv with constant denominator in terms of other ids. After extracting these divisions, local ids in other
with the same division representation as some other local id in any set are considered duplicate and are merged.
It is possible that division representation for some local id cannot be obtained, and thus these local ids are not considered for detecting duplicates.
Definition at line 1356 of file IntegerRelation.cpp.
References eliminateRedundantLocalVar(), getNumLocalVars(), and mlir::presburger::mergeLocalVars().
Referenced by mlir::affine::MemRefAccess::getAccessRelation(), getSetDifference(), intersect(), intersectDomain(), intersectRange(), mergeAndAlignVars(), and mergeAndCompose().
|
protected |
Normalized each constraints by the GCD of its coefficients.
Definition at line 585 of file IntegerRelation.cpp.
References equalities, getNumEqualities(), getNumInequalities(), inequalities, and mlir::presburger::IntMatrix::normalizeRow().
Referenced by projectOut(), removeTrivialRedundancy(), and simplify().
void IntegerRelation::print | ( | raw_ostream & | os | ) | const |
Definition at line 2592 of file IntegerRelation.cpp.
References atEq(), atIneq(), getNumCols(), getNumEqualities(), getNumInequalities(), hasConsistentState(), and printSpace().
Referenced by dump(), and mlir::presburger::PresburgerRelation::print().
|
protectedvirtual |
Prints the number of constraints, dimensions, symbols and locals in the IntegerRelation.
Reimplemented in mlir::FlatLinearValueConstraints, and mlir::FlatLinearConstraints.
Definition at line 2503 of file IntegerRelation.cpp.
References getNumConstraints(), mlir::presburger::PresburgerSpace::print(), and space.
Referenced by print().
|
inline |
Definition at line 479 of file IntegerRelation.h.
References projectOut().
Referenced by projectOut().
void IntegerRelation::projectOut | ( | unsigned | pos, |
unsigned | num | ||
) |
Projects out (aka eliminates) num
variables starting at position pos
.
The resulting constraint system is the shadow along the dimensions that still exist. This method may not always be integer exact.
Definition at line 2096 of file IntegerRelation.cpp.
References fourierMotzkinEliminate(), gaussianEliminateVars(), gcdTightenInequalities(), getBestVarToEliminate(), getNumCols(), getNumVars(), and normalizeConstraintsByGCD().
Referenced by mlir::ValueBoundsConstraintSet::computeBound(), computeConstantLowerOrUpperBound(), computeDirectionVector(), and mlir::ValueBoundsConstraintSet::projectOut().
|
protected |
Checks for identical inequalities and eliminates redundant inequalities.
Returns whether the constraint system was modified.
Definition at line 2360 of file IntegerRelation.cpp.
References addEquality(), atIneq(), cols, getEmpty(), getInequality(), getNumCols(), getNumInequalities(), getSpace(), inequalities, removeInequality(), and mlir::presburger::Matrix< T >::swapRows().
Referenced by simplify().
void IntegerRelation::removeDuplicateDivs | ( | ) |
Definition at line 1391 of file IntegerRelation.cpp.
References eliminateRedundantLocalVar(), getLocalReprs(), and mlir::presburger::DivisionRepr::removeDuplicateDivs().
Referenced by getSetDifference().
void IntegerRelation::removeEquality | ( | unsigned | pos | ) |
Definition at line 442 of file IntegerRelation.cpp.
References equalities, and mlir::presburger::Matrix< T >::removeRow().
Referenced by gaussianEliminateVars(), removeConstraintsInvolvingVarRange(), removeIndependentConstraints(), removeRedundantLocalVars(), and removeTrivialEqualities().
void IntegerRelation::removeEqualityRange | ( | unsigned | start, |
unsigned | end | ||
) |
Remove the (in)equalities at positions [start, end).
Definition at line 450 of file IntegerRelation.cpp.
References equalities, and mlir::presburger::Matrix< T >::removeRows().
Referenced by gaussianEliminate(), and truncate().
void IntegerRelation::removeIndependentConstraints | ( | unsigned | pos, |
unsigned | num | ||
) |
Removes constraints that are independent of (i.e., do not have a coefficient) variables in the range [pos, pos + num).
Definition at line 2330 of file IntegerRelation.cpp.
References getIndependentConstraints(), getNumVars(), removeEquality(), and removeInequality().
Referenced by createSeparationCondition().
void IntegerRelation::removeInequality | ( | unsigned | pos | ) |
Definition at line 446 of file IntegerRelation.cpp.
References inequalities, and mlir::presburger::Matrix< T >::removeRow().
Referenced by mlir::ValueBoundsConstraintSet::comparePos(), mlir::presburger::MultiAffineFunction::getLexSet(), removeConstraintsInvolvingVarRange(), removeDuplicateConstraints(), and removeIndependentConstraints().
void IntegerRelation::removeInequalityRange | ( | unsigned | start, |
unsigned | end | ||
) |
Definition at line 456 of file IntegerRelation.cpp.
References inequalities, and mlir::presburger::Matrix< T >::removeRows().
Referenced by truncate().
void IntegerRelation::removeRedundantConstraints | ( | ) |
Removes redundant constraints using Simplex.
Although the algorithm can theoretically take exponential time in the worst case (rare), it is known to perform much better in the average case. If V is the number of vertices in the polytope and C is the number of constraints, the algorithm takes O(VC) time.
Definition at line 1210 of file IntegerRelation.cpp.
References mlir::presburger::Simplex::detectRedundant(), gcdTightenInequalities(), and getNumInequalities().
void IntegerRelation::removeRedundantInequalities | ( | ) |
A more expensive check than removeTrivialRedundancy
to detect redundant inequalities.
Definition at line 1178 of file IntegerRelation.cpp.
References getNumInequalities().
|
protected |
Removes local variables using equalities.
Each equality is checked if it can be reduced to the form: e = affine-expr
, where e
is a local variable and affine-expr
is an affine expression not containing e
. If an equality satisfies this form, the local variable is replaced in each constraint and then removed. The equality used to replace this local variable is also removed.
Definition at line 1420 of file IntegerRelation.cpp.
References mlir::presburger::abs(), atEq(), eliminateFromConstraint(), equalities, getNumDimAndSymbolVars(), getNumEqualities(), getNumInequalities(), getNumVars(), mlir::presburger::IntMatrix::normalizeRow(), removeEquality(), and removeVar().
Referenced by mlir::affine::FlatAffineRelation::compose(), isEmpty(), and mergeAndCompose().
void IntegerRelation::removeTrivialEqualities | ( | ) |
Definition at line 2508 of file IntegerRelation.cpp.
References getEquality(), getNumEqualities(), rangeIsZero(), and removeEquality().
Referenced by isFullDim().
void IntegerRelation::removeTrivialRedundancy | ( | ) |
Removes duplicate constraints, trivially true constraints, and constraints that can be detected as redundant as a result of differing only in their constant term part.
A constraint of the form <non-negative constant> >= 0 is considered trivially true. This method is a linear time method on the constraints, does a single scan, and updates in place. It also normalizes constraints by their GCD and performs GCD tightening on inequalities.
A constraint of the form <non-negative constant> >= 0 is considered trivially true.
Definition at line 1828 of file IntegerRelation.cpp.
References gcdTightenInequalities(), and normalizeConstraintsByGCD().
Referenced by createSeparationCondition(), mlir::affine::simplifyIntegerSet(), and unionBoundingBox().
void IntegerRelation::removeVar | ( | unsigned | pos | ) |
Removes the specified variable from the system.
Definition at line 386 of file IntegerRelation.cpp.
References removeVarRange().
void IntegerRelation::removeVar | ( | VarKind | kind, |
unsigned | pos | ||
) |
Removes variables of the specified kind with the specified pos (or within the specified range) from the system.
The specified location is relative to the first variable of the specified kind.
Definition at line 382 of file IntegerRelation.cpp.
References removeVarRange().
Referenced by createSeparationCondition(), eliminateRedundantLocalVar(), and removeRedundantLocalVars().
|
protected |
Removes variables in the column range [varStart, varLimit), and copies any remaining valid data into place, updates member variables, and resizes arrays as needed.
Definition at line 404 of file IntegerRelation.cpp.
References mlir::presburger::Domain, getNumVarKind(), getNumVars(), getVarKindOffset(), mlir::presburger::Local, min(), mlir::presburger::Range, removeVarRange(), and mlir::presburger::Symbol.
|
virtual |
Reimplemented in mlir::affine::FlatAffineRelation, and mlir::FlatLinearValueConstraints.
Definition at line 388 of file IntegerRelation.cpp.
References equalities, getNumVarKind(), getVarKindOffset(), inequalities, mlir::presburger::Matrix< T >::removeColumns(), mlir::presburger::PresburgerSpace::removeVarRange(), and space.
Referenced by gaussianEliminateVars(), mergeAndCompose(), removeVar(), removeVarRange(), setAndEliminate(), and truncateVarKind().
|
inline |
Definition at line 138 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::resetIds(), and space.
Referenced by mlir::affine::MemRefAccess::getAccessRelation().
void IntegerRelation::setAndEliminate | ( | unsigned | pos, |
ArrayRef< DynamicAPInt > | values | ||
) |
Sets the values.size()
variables starting at po
s to the specified values and removes them.
Definition at line 546 of file IntegerRelation.cpp.
References mlir::presburger::Matrix< T >::addToColumn(), equalities, getNumCols(), getNumVars(), inequalities, and removeVarRange().
Referenced by constantFoldVar(), and setAndEliminate().
|
inline |
Definition at line 366 of file IntegerRelation.h.
References mlir::presburger::getDynamicAPIntVec(), and setAndEliminate().
|
inline |
Changes the partition between dimensions and symbols.
Depending on the new symbol count, either a chunk of dimensional variables immediately before the split become symbols, or some of the symbols immediately after the split become dimensions.
Definition at line 649 of file IntegerRelation.h.
References mlir::presburger::PresburgerSpace::setVarSymbolSeparation(), and space.
Referenced by createFullTiles(), createSeparationCondition(), and turnSymbolIntoDim().
void IntegerRelation::setId | ( | VarKind | kind, |
unsigned | i, | ||
Identifier | id | ||
) |
Set the identifier for the ith variable of the specified kind of the IntegerRelation's PresburgerSpace.
The index is relative to the kind of the variable.
Definition at line 71 of file IntegerRelation.cpp.
References mlir::presburger::PresburgerSpace::getNumVarKind(), mlir::presburger::PresburgerSpace::isUsingIds(), mlir::presburger::Local, mlir::presburger::PresburgerSpace::setId(), and space.
Referenced by mlir::affine::MemRefAccess::getAccessRelation(), and mlir::affine::getRelationFromMap().
void IntegerRelation::setSpace | ( | const PresburgerSpace & | oSpace | ) |
Set the space to oSpace
, which should have the same number of ids as the current space.
Definition at line 58 of file IntegerRelation.cpp.
References mlir::presburger::PresburgerSpace::getNumVars(), and space.
void IntegerRelation::setSpaceExceptLocals | ( | const PresburgerSpace & | oSpace | ) |
Set the space to oSpace
, which should not have any local ids.
oSpace
can have fewer ids than the current space; in that case, the the extra ids in this
that are not accounted for by oSpace
will be considered as local ids. oSpace
should not have more ids than the current space; this will result in an assert failure.
Definition at line 63 of file IntegerRelation.cpp.
References mlir::presburger::PresburgerSpace::getNumLocalVars(), getNumVars(), mlir::presburger::PresburgerSpace::getNumVars(), mlir::presburger::PresburgerSpace::insertVar(), mlir::presburger::Local, and space.
void IntegerRelation::simplify | ( | ) |
Simplify the constraint system by removing canonicalizing constraints and removing redundant constraints.
Definition at line 1400 of file IntegerRelation.cpp.
References gaussianEliminate(), isObviouslyEmpty(), normalizeConstraintsByGCD(), and removeDuplicateConstraints().
Referenced by mlir::presburger::PresburgerRelation::simplify().
PresburgerRelation IntegerRelation::subtract | ( | const PresburgerRelation & | set | ) | const |
Return the set difference of this set and the given set, i.e., return this \ set
.
Definition at line 350 of file IntegerRelation.cpp.
References mlir::presburger::PresburgerRelation::subtract().
Referenced by mlir::presburger::IntegerPolyhedron::subtract().
|
virtual |
Swap the posA^th variable with the posB^th variable.
Definition at line 462 of file IntegerRelation.cpp.
References equalities, getNumVars(), mlir::presburger::PresburgerSpace::getVarKindAt(), getVarKindOffset(), inequalities, space, mlir::presburger::Matrix< T >::swapColumns(), and mlir::presburger::PresburgerSpace::swapVar().
Referenced by mlir::affine::MemRefAccess::getAccessRelation(), mergeAndAlignSymbols(), mergeAndAlignVars(), mlir::FlatLinearValueConstraints::mergeSymbolVars(), and turnSymbolIntoDim().
void IntegerRelation::truncate | ( | const CountsSnapshot & | counts | ) |
Definition at line 214 of file IntegerRelation.cpp.
References mlir::presburger::Domain, mlir::presburger::IntegerRelation::CountsSnapshot::getNumEqs(), getNumEqualities(), mlir::presburger::IntegerRelation::CountsSnapshot::getNumIneqs(), getNumInequalities(), mlir::presburger::Local, mlir::presburger::Range, removeEqualityRange(), removeInequalityRange(), mlir::presburger::Symbol, and truncateVarKind().
Referenced by mlir::presburger::SymbolicLexSimplex::computeSymbolicIntegerLexMin(), and getSetDifference().
|
protected |
Truncate the vars to the number in the space of the specified CountsSnapshot.
Definition at line 209 of file IntegerRelation.cpp.
References mlir::presburger::PresburgerSpace::getNumVarKind(), mlir::presburger::IntegerRelation::CountsSnapshot::getSpace(), and truncateVarKind().
|
protected |
Truncate the vars of the specified kind to the specified number by dropping some vars at the end.
num
must be less than the current number.
Definition at line 203 of file IntegerRelation.cpp.
References getNumVarKind(), and removeVarRange().
Referenced by truncate(), and truncateVarKind().
LogicalResult IntegerRelation::unionBoundingBox | ( | const IntegerRelation & | other | ) |
Updates the constraints to be the smallest bounding (enclosing) box that contains the points of this
set and that of other
, with the symbols being treated specially.
For each of the dimensions, the min of the lower bounds (symbolic) and the max of the upper bounds (symbolic) is computed to determine such a bounding box. other
is expected to have the same dimensional variables as this constraint system (in the same order).
E.g.: 1) this = {0 <= d0 <= 127}, other = {16 <= d0 <= 192}, output = {0 <= d0 <= 192} 2) this = {s0 + 5 <= d0 <= s0 + 20}, other = {s0 + 1 <= d0 <= s0 + 9}, output = {s0 + 1 <= d0 <= s0 + 20} 3) this = {0 <= d0 <= 5, 1 <= d1 <= 9} other = {2 <= d0 <= 6, 5 <= d1 <= 15}, output = {0 <= d0 <= 6, 1 <= d1 <= 15}
Definition at line 2181 of file IntegerRelation.cpp.
References addInequality(), append(), clearConstraints(), copy(), getCommonConstraints(), getConstantBound(), getConstantBoundOnDimSize(), getNumCols(), getNumDimVars(), getNumLocalVars(), getNumSymbolVars(), mlir::presburger::PresburgerSpace::getRelationSpace(), getSpace(), mlir::presburger::PresburgerSpace::isEqual(), mlir::presburger::LB, max(), min(), removeTrivialRedundancy(), space, and mlir::presburger::UB.
|
protected |
Coefficients of affine equalities (in == 0 form).
Definition at line 874 of file IntegerRelation.h.
Referenced by addBound(), addEquality(), append(), atEq(), atEq64(), clearConstraints(), convertVarKind(), eliminateRedundantLocalVar(), gaussianEliminate(), gaussianEliminateVars(), getEquality(), getEquality64(), getNumEqualities(), getNumReservedEqualities(), hasConsistentState(), insertVar(), normalizeConstraintsByGCD(), removeEquality(), removeEqualityRange(), removeRedundantLocalVars(), removeVarRange(), setAndEliminate(), and swapVar().
|
protected |
Coefficients of affine inequalities (in >= 0 form).
Definition at line 877 of file IntegerRelation.h.
Referenced by addBound(), addInequality(), append(), atIneq(), atIneq64(), clearConstraints(), convertVarKind(), eliminateRedundantLocalVar(), gaussianEliminate(), gaussianEliminateVars(), gcdTightenInequalities(), getInequalities(), getInequality(), getInequality64(), getNumInequalities(), getNumReservedInequalities(), hasConsistentState(), insertVar(), mlir::presburger::IntegerPolyhedron::IntegerPolyhedron(), normalizeConstraintsByGCD(), removeDuplicateConstraints(), removeInequality(), removeInequalityRange(), removeVarRange(), setAndEliminate(), and swapVar().
|
staticconstexprprotected |
A parameter that controls detection of an unrealistic number of constraints.
If the number of constraints is this many times the number of variables, we consider such a system out of line with the intended use case of IntegerRelation.
Definition at line 869 of file IntegerRelation.h.
Referenced by isEmpty().
|
protected |
Definition at line 871 of file IntegerRelation.h.
Referenced by append(), mlir::FlatLinearValueConstraints::computeAlignedMap(), computeReprWithOnlyDivLocals(), convertVarKind(), getEmpty(), getIds(), mlir::FlatLinearValueConstraints::getMaybeValues(), getNumCols(), getNumDimAndSymbolVars(), getNumDimVars(), getNumDomainVars(), getNumLocalVars(), getNumRangeVars(), getNumSymbolVars(), getNumVarKind(), getNumVars(), getSpace(), getSpaceWithoutLocals(), getUniverse(), mlir::presburger::IntegerPolyhedron::getUniverse(), mlir::FlatLinearValueConstraints::getValue(), getVarKindAt(), getVarKindEnd(), getVarKindOffset(), getVarKindOverlap(), mlir::FlatLinearValueConstraints::hasValue(), insertVar(), mlir::presburger::IntegerPolyhedron::IntegerPolyhedron(), IntegerRelation(), isEqual(), isObviouslyEqual(), isSubsetOf(), mergeAndAlignSymbols(), mergeAndCompose(), printSpace(), removeVarRange(), resetIds(), setDimSymbolSeparation(), setId(), setSpace(), setSpaceExceptLocals(), mlir::FlatLinearValueConstraints::setValue(), swapVar(), and unionBoundingBox().