My Project
Public Member Functions | Private Member Functions | Private Attributes
tropicalStrategy Class Reference

#include <tropicalStrategy.h>

Public Member Functions

 tropicalStrategy (const ideal I, const ring r, const bool completelyHomogeneous=true, const bool completeSpace=true)
 Constructor for the trivial valuation case. More...
 
 tropicalStrategy (const ideal J, const number p, const ring s)
 Constructor for the non-trivial valuation case p is the uniformizing parameter of the valuation. More...
 
 tropicalStrategy (const tropicalStrategy &currentStrategy)
 copy constructor More...
 
 ~tropicalStrategy ()
 destructor More...
 
tropicalStrategyoperator= (const tropicalStrategy &currentStrategy)
 assignment operator More...
 
bool isValuationTrivial () const
 
bool isValuationNonTrivial () const
 
ring getOriginalRing () const
 returns the polynomial ring over the field with valuation More...
 
ideal getOriginalIdeal () const
 returns the input ideal over the field with valuation More...
 
ring getStartingRing () const
 returns the polynomial ring over the valuation ring More...
 
ideal getStartingIdeal () const
 returns the input ideal More...
 
int getExpectedAmbientDimension () const
 
int getExpectedDimension () const
 returns the expected Dimension of the polyhedral output More...
 
number getUniformizingParameter () const
 returns the uniformizing parameter in the valuation ring More...
 
ring getShortcutRing () const
 
gfan::ZCone getHomogeneitySpace () const
 returns the homogeneity space of the preimage ideal More...
 
bool homogeneitySpaceContains (const gfan::ZVector &v) const
 returns true, if v is contained in the homogeneity space; false otherwise More...
 
bool restrictToLowerHalfSpace () const
 returns true, if valuation non-trivial, false otherwise More...
 
gfan::ZVector adjustWeightForHomogeneity (gfan::ZVector w) const
 Given weight w, returns a strictly positive weight u such that an ideal satisfying the valuation-sepcific homogeneity conditions is weighted homogeneous with respect to w if and only if it is homogeneous with respect to u. More...
 
gfan::ZVector adjustWeightUnderHomogeneity (gfan::ZVector v, gfan::ZVector w) const
 Given strictly positive weight w and weight v, returns a strictly positive weight u such that on an ideal that is weighted homogeneous with respect to w the weights u and v coincide. More...
 
gfan::ZVector negateWeight (const gfan::ZVector &w) const
 
ring getShortcutRingPrependingWeight (const ring r, const gfan::ZVector &w) const
 If valuation trivial, returns a copy of r with a positive weight prepended, such that any ideal homogeneous with respect to w is homogeneous with respect to that weight. More...
 
bool reduce (ideal I, const ring r) const
 reduces the generators of an ideal I so that the inequalities and equations of the Groebner cone can be read off. More...
 
void pReduce (ideal I, const ring r) const
 
std::pair< poly, int > checkInitialIdealForMonomial (const ideal I, const ring r, const gfan::ZVector &w=0) const
 If given w, assuming w is in the Groebner cone of the ordering on r and I is a standard basis with respect to that ordering, checks whether the initial ideal of I with respect to w contains a monomial. More...
 
ideal computeStdOfInitialIdeal (const ideal inI, const ring r) const
 given generators of the initial ideal, computes its standard basis More...
 
ideal computeWitness (const ideal inJ, const ideal inI, const ideal I, const ring r) const
 suppose w a weight in maximal groebner cone of > suppose I (initially) reduced standard basis w.r.t. More...
 
ideal computeLift (const ideal inJs, const ring s, const ideal inIr, const ideal Ir, const ring r) const
 
std::pair< ideal, ring > computeFlip (const ideal Ir, const ring r, const gfan::ZVector &interiorPoint, const gfan::ZVector &facetNormal) const
 given an interior point of a groebner cone computes the groebner cone adjacent to it More...
 

Private Member Functions

ring copyAndChangeCoefficientRing (const ring r) const
 
ring copyAndChangeOrderingWP (const ring r, const gfan::ZVector &w, const gfan::ZVector &v) const
 
ring copyAndChangeOrderingLS (const ring r, const gfan::ZVector &w, const gfan::ZVector &v) const
 
bool checkForUniformizingBinomial (const ideal I, const ring r) const
 if valuation non-trivial, checks whether the generating system contains p-t otherwise returns true More...
 
bool checkForUniformizingParameter (const ideal inI, const ring r) const
 if valuation non-trivial, checks whether the genearting system contains p otherwise returns true More...
 
int findPositionOfUniformizingBinomial (const ideal I, const ring r) const
 
void putUniformizingBinomialInFront (ideal I, const ring r, const number q) const
 

Private Attributes

ring originalRing
 polynomial ring over a field with valuation More...
 
ideal originalIdeal
 input ideal, assumed to be a homogeneous prime ideal More...
 
int expectedDimension
 the expected Dimension of the polyhedral output, i.e. More...
 
gfan::ZCone linealitySpace
 the homogeneity space of the Grobner fan More...
 
ring startingRing
 polynomial ring over the valuation ring extended by one extra variable t More...
 
ideal startingIdeal
 preimage of the input ideal under the map that sends t to the uniformizing parameter More...
 
number uniformizingParameter
 uniformizing parameter in the valuation ring More...
 
ring shortcutRing
 polynomial ring over the residue field More...
 
bool onlyLowerHalfSpace
 true if valuation non-trivial, false otherwise More...
 
gfan::ZVector(* weightAdjustingAlgorithm1 )(const gfan::ZVector &w)
 A function such that: Given weight w, returns a strictly positive weight u such that an ideal satisfying the valuation-sepcific homogeneity conditions is weighted homogeneous with respect to w if and only if it is homogeneous with respect to u. More...
 
gfan::ZVector(* weightAdjustingAlgorithm2 )(const gfan::ZVector &v, const gfan::ZVector &w)
 A function such that: Given strictly positive weight w and weight v, returns a strictly positive weight u such that on an ideal that is weighted homogeneous with respect to w the weights u and v coincide. More...
 
bool(* extraReductionAlgorithm )(ideal I, ring r, number p)
 A function that reduces the generators of an ideal I so that the inequalities and equations of the Groebner cone can be read off. More...
 

Detailed Description

Definition at line 36 of file tropicalStrategy.h.

Constructor & Destructor Documentation

◆ tropicalStrategy() [1/3]

tropicalStrategy::tropicalStrategy ( const ideal  I,
const ring  r,
const bool  completelyHomogeneous = true,
const bool  completeSpace = true 
)

Constructor for the trivial valuation case.

Initializes all relevant structures and information for the trivial valuation case, i.e.

computing a tropical variety without any valuation.

Definition at line 138 of file tropicalStrategy.cc.

140 :
149 onlyLowerHalfSpace(false),
153{
155 if (!completelyHomogeneous)
156 {
159 }
160 if (!completeSpace)
161 onlyLowerHalfSpace = true;
162}
gfan::ZVector nonvalued_adjustWeightForHomogeneity(const gfan::ZVector &w)
gfan::ZVector nonvalued_adjustWeightUnderHomogeneity(const gfan::ZVector &e, const gfan::ZVector &)
gfan::ZVector valued_adjustWeightForHomogeneity(const gfan::ZVector &w)
gfan::ZVector valued_adjustWeightUnderHomogeneity(const gfan::ZVector &e, const gfan::ZVector &w)
int expectedDimension
the expected Dimension of the polyhedral output, i.e.
bool onlyLowerHalfSpace
true if valuation non-trivial, false otherwise
gfan::ZCone linealitySpace
the homogeneity space of the Grobner fan
ring startingRing
polynomial ring over the valuation ring extended by one extra variable t
ideal originalIdeal
input ideal, assumed to be a homogeneous prime ideal
gfan::ZVector(* weightAdjustingAlgorithm1)(const gfan::ZVector &w)
A function such that: Given weight w, returns a strictly positive weight u such that an ideal satisfy...
ring shortcutRing
polynomial ring over the residue field
bool(* extraReductionAlgorithm)(ideal I, ring r, number p)
A function that reduces the generators of an ideal I so that the inequalities and equations of the Gr...
number uniformizingParameter
uniformizing parameter in the valuation ring
ideal startingIdeal
preimage of the input ideal under the map that sends t to the uniformizing parameter
gfan::ZVector(* weightAdjustingAlgorithm2)(const gfan::ZVector &v, const gfan::ZVector &w)
A function such that: Given strictly positive weight w and weight v, returns a strictly positive weig...
ring originalRing
polynomial ring over a field with valuation
ideal id_Copy(ideal h1, const ring r)
copy an ideal
#define assume(x)
Definition: mod2.h:389
#define NULL
Definition: omList.c:12
ring rCopy(ring r)
Definition: ring.cc:1731
static BOOLEAN rField_is_Z(const ring r)
Definition: ring.h:509
static BOOLEAN rField_is_Zp(const ring r)
Definition: ring.h:500
static BOOLEAN rField_is_Q(const ring r)
Definition: ring.h:506
static bool noExtraReduction(ideal I, ring r, number)
int dim(ideal I, ring r)
gfan::ZCone homogeneitySpace(ideal I, ring r)
Definition: tropical.cc:19

◆ tropicalStrategy() [2/3]

tropicalStrategy::tropicalStrategy ( const ideal  J,
const number  p,
const ring  s 
)

Constructor for the non-trivial valuation case p is the uniformizing parameter of the valuation.

Definition at line 277 of file tropicalStrategy.cc.

277 :
281 linealitySpace(gfan::ZCone()), // to come, see below
282 startingRing(NULL), // to come, see below
283 startingIdeal(NULL), // to come, see below
284 uniformizingParameter(NULL), // to come, see below
285 shortcutRing(NULL), // to come, see below
286 onlyLowerHalfSpace(true),
290{
291 /* assume that the ground field of the originalRing is Q */
293
294 /* replace Q with Z for the startingRing
295 * and add an extra variable for tracking the uniformizing parameter */
297
298 /* map the uniformizing parameter into the new coefficient domain */
301
302 /* map the input ideal into the new polynomial ring */
305
307
308 /* construct the shorcut ring */
309 shortcutRing = rCopy0(startingRing,FALSE); // do not copy q-ideal
314}
#define FALSE
Definition: auxiliary.h:96
bool reduce(ideal I, const ring r) const
reduces the generators of an ideal I so that the inequalities and equations of the Groebner cone can ...
static FORCE_INLINE long n_Int(number &n, const coeffs r)
conversion of n to an int; 0 if not possible in Z/pZ: the representing int lying in (-p/2 ....
Definition: coeffs.h:544
@ n_Zp
\F{p < 2^31}
Definition: coeffs.h:29
static FORCE_INLINE nMapFunc n_SetMap(const coeffs src, const coeffs dst)
set the mapping function pointers for translating numbers from src to dst
Definition: coeffs.h:697
coeffs nInitChar(n_coeffType t, void *parameter)
one-time initialisations for new coeffs in case of an error return NULL
Definition: numbers.cc:413
number(* nMapFunc)(number a, const coeffs src, const coeffs dst)
maps "a", which lives in src, into dst
Definition: coeffs.h:73
void nKillChar(coeffs r)
undo all initialisations
Definition: numbers.cc:568
const CanonicalForm int s
Definition: facAbsFact.cc:51
bool ppreduceInitially(poly *hStar, const poly g, const ring r)
reduces h initially with respect to g, returns false if h was initially reduced in the first place,...
int IsPrime(int p)
Definition: prime.cc:61
BOOLEAN rComplete(ring r, int force)
this needs to be called whenever a new ring is created: new fields in ring are created (like VarOffse...
Definition: ring.cc:3450
ring rCopy0(const ring r, BOOLEAN copy_qideal, BOOLEAN copy_ordering)
Definition: ring.cc:1421
#define rTest(r)
Definition: ring.h:782
static ideal constructStartingIdeal(ideal originalIdeal, ring originalRing, number uniformizingParameter, ring startingRing)
static ring constructStartingRing(ring r)
Given a polynomial ring r over the rational numbers and a weighted ordering, returns a polynomial rin...

◆ tropicalStrategy() [3/3]

tropicalStrategy::tropicalStrategy ( const tropicalStrategy currentStrategy)

copy constructor

Definition at line 316 of file tropicalStrategy.cc.

316 :
317 originalRing(rCopy(currentStrategy.getOriginalRing())),
318 originalIdeal(id_Copy(currentStrategy.getOriginalIdeal(),currentStrategy.getOriginalRing())),
319 expectedDimension(currentStrategy.getExpectedDimension()),
320 linealitySpace(currentStrategy.getHomogeneitySpace()),
321 startingRing(rCopy(currentStrategy.getStartingRing())),
322 startingIdeal(id_Copy(currentStrategy.getStartingIdeal(),currentStrategy.getStartingRing())),
329{
334 if (currentStrategy.getUniformizingParameter())
335 {
338 }
339 if (currentStrategy.getShortcutRing())
340 {
341 shortcutRing = rCopy(currentStrategy.getShortcutRing());
343 }
344}
ideal getOriginalIdeal() const
returns the input ideal over the field with valuation
gfan::ZCone getHomogeneitySpace() const
returns the homogeneity space of the preimage ideal
int getExpectedDimension() const
returns the expected Dimension of the polyhedral output
ring getStartingRing() const
returns the polynomial ring over the valuation ring
ideal getStartingIdeal() const
returns the input ideal
bool restrictToLowerHalfSpace() const
returns true, if valuation non-trivial, false otherwise
ring getOriginalRing() const
returns the polynomial ring over the field with valuation
number getUniformizingParameter() const
returns the uniformizing parameter in the valuation ring
ring getShortcutRing() const
static FORCE_INLINE number n_Copy(number n, const coeffs r)
return a copy of 'n'
Definition: coeffs.h:448
#define n_Test(a, r)
BOOLEAN n_Test(number a, const coeffs r)
Definition: coeffs.h:709
#define id_Test(A, lR)
Definition: simpleideals.h:87

◆ ~tropicalStrategy()

tropicalStrategy::~tropicalStrategy ( )

destructor

Definition at line 346 of file tropicalStrategy.cc.

347{
354
361}
static FORCE_INLINE void n_Delete(number *p, const coeffs r)
delete 'p'
Definition: coeffs.h:452
void rDelete(ring r)
unconditionally deletes fields in r
Definition: ring.cc:450
void id_Delete(ideal *h, ring r)
deletes an ideal/module/matrix

Member Function Documentation

◆ adjustWeightForHomogeneity()

gfan::ZVector tropicalStrategy::adjustWeightForHomogeneity ( gfan::ZVector  w) const
inline

Given weight w, returns a strictly positive weight u such that an ideal satisfying the valuation-sepcific homogeneity conditions is weighted homogeneous with respect to w if and only if it is homogeneous with respect to u.

Definition at line 248 of file tropicalStrategy.h.

249 {
250 return this->weightAdjustingAlgorithm1(w);
251 }
const CanonicalForm & w
Definition: facAbsFact.cc:51

◆ adjustWeightUnderHomogeneity()

gfan::ZVector tropicalStrategy::adjustWeightUnderHomogeneity ( gfan::ZVector  v,
gfan::ZVector  w 
) const
inline

Given strictly positive weight w and weight v, returns a strictly positive weight u such that on an ideal that is weighted homogeneous with respect to w the weights u and v coincide.

Definition at line 258 of file tropicalStrategy.h.

259 {
260 return this->weightAdjustingAlgorithm2(v,w);
261 }
const Variable & v
< [in] a sqrfree bivariate poly
Definition: facBivar.h:39

◆ checkForUniformizingBinomial()

bool tropicalStrategy::checkForUniformizingBinomial ( const ideal  I,
const ring  r 
) const
private

if valuation non-trivial, checks whether the generating system contains p-t otherwise returns true

Definition at line 814 of file tropicalStrategy.cc.

815{
816 // if the valuation is trivial,
817 // then there is no special condition the first generator has to fullfill
818 if (isValuationTrivial())
819 return true;
820
821 // if the valuation is non-trivial then checks if the first generator is p-t
822 nMapFunc identity = n_SetMap(startingRing->cf,r->cf);
823 poly p = p_One(r);
824 p_SetCoeff(p,identity(uniformizingParameter,startingRing->cf,r->cf),r);
825 poly t = p_One(r);
826 p_SetExp(t,1,1,r);
827 p_Setm(t,r);
828 poly pt = p_Add_q(p,p_Neg(t,r),r);
829
830 for (int i=0; i<IDELEMS(I); i++)
831 {
832 if (p_EqualPolys(I->m[i],pt,r))
833 {
834 p_Delete(&pt,r);
835 return true;
836 }
837 }
838 p_Delete(&pt,r);
839 return false;
840}
int i
Definition: cfEzgcd.cc:132
int p
Definition: cfModGcd.cc:4078
bool isValuationTrivial() const
poly p_One(const ring r)
Definition: p_polys.cc:1313
BOOLEAN p_EqualPolys(poly p1, poly p2, const ring r)
Definition: p_polys.cc:4508
static poly p_Neg(poly p, const ring r)
Definition: p_polys.h:1105
static poly p_Add_q(poly p, poly q, const ring r)
Definition: p_polys.h:934
static unsigned long p_SetExp(poly p, const unsigned long e, const unsigned long iBitmask, const int VarOffset)
set a single variable exponent @Note: VarOffset encodes the position in p->exp
Definition: p_polys.h:486
static void p_Setm(poly p, const ring r)
Definition: p_polys.h:231
static number p_SetCoeff(poly p, number n, ring r)
Definition: p_polys.h:410
static void p_Delete(poly *p, const ring r)
Definition: p_polys.h:899
#define IDELEMS(i)
Definition: simpleideals.h:23

◆ checkForUniformizingParameter()

bool tropicalStrategy::checkForUniformizingParameter ( const ideal  inI,
const ring  r 
) const
private

if valuation non-trivial, checks whether the genearting system contains p otherwise returns true

Definition at line 867 of file tropicalStrategy.cc.

868{
869 // if the valuation is trivial,
870 // then there is no special condition the first generator has to fullfill
871 if (isValuationTrivial())
872 return true;
873
874 // if the valuation is non-trivial then checks if the first generator is p
875 if (inI->m[0]==NULL)
876 return false;
877 nMapFunc identity = n_SetMap(startingRing->cf,r->cf);
878 poly p = p_One(r);
879 p_SetCoeff(p,identity(uniformizingParameter,startingRing->cf,r->cf),r);
880
881 for (int i=0; i<IDELEMS(inI); i++)
882 {
883 if (p_EqualPolys(inI->m[i],p,r))
884 {
885 p_Delete(&p,r);
886 return true;
887 }
888 }
889 p_Delete(&p,r);
890 return false;
891}

◆ checkInitialIdealForMonomial()

std::pair< poly, int > tropicalStrategy::checkInitialIdealForMonomial ( const ideal  I,
const ring  r,
const gfan::ZVector &  w = 0 
) const

If given w, assuming w is in the Groebner cone of the ordering on r and I is a standard basis with respect to that ordering, checks whether the initial ideal of I with respect to w contains a monomial.

If no w is given, assuming that I is already an initial form of some ideal, checks whether I contains a monomial. In both cases returns a monomial, if it contains one, returns NULL otherwise.

Definition at line 496 of file tropicalStrategy.cc.

497{
498 // quick check whether I already contains an ideal
499 int k = IDELEMS(I);
500 for (int i=0; i<k; i++)
501 {
502 poly g = I->m[i];
503 if (g!=NULL
504 && pNext(g)==NULL
505 && (isValuationTrivial() || n_IsOne(p_GetCoeff(g,r),r->cf)))
506 return std::pair<poly,int>(g,i);
507 }
508
509 ring rShortcut;
510 ideal inIShortcut;
511 if (w.size()>0)
512 {
513 // if needed, prepend extra weight for homogeneity
514 // switch to residue field if valuation is non trivial
515 rShortcut = getShortcutRingPrependingWeight(r,w);
516
517 // compute the initial ideal and map it into the constructed ring
518 // if switched to residue field, remove possibly 0 elements
519 ideal inI = initial(I,r,w);
520 inIShortcut = idInit(k);
521 nMapFunc intoShortcut = n_SetMap(r->cf,rShortcut->cf);
522 for (int i=0; i<k; i++)
523 inIShortcut->m[i] = p_PermPoly(inI->m[i],NULL,r,rShortcut,intoShortcut,NULL,0);
525 idSkipZeroes(inIShortcut);
526 id_Delete(&inI,r);
527 }
528 else
529 {
530 rShortcut = r;
531 inIShortcut = I;
532 }
533
534 gfan::ZCone C0 = homogeneitySpace(inIShortcut,rShortcut);
535 gfan::ZCone pos = gfan::ZCone::positiveOrthant(C0.ambientDimension());
536 gfan::ZCone C0pos = intersection(C0,pos);
537 C0pos.canonicalize();
538 gfan::ZVector wpos = C0pos.getRelativeInteriorPoint();
540
541 // check initial ideal for monomial and
542 // if it exsists, return a copy of the monomial in the input ring
543 poly p = searchForMonomialViaStepwiseSaturation(inIShortcut,rShortcut,wpos);
544 poly monomial = NULL;
545 if (p!=NULL)
546 {
547 monomial=p_One(r);
548 for (int i=1; i<=rVar(r); i++)
549 p_SetExp(monomial,i,p_GetExp(p,i,rShortcut),r);
550 p_Setm(monomial,r);
551 p_Delete(&p,rShortcut);
552 }
553
554
555 if (w.size()>0)
556 {
557 // if needed, cleanup
558 id_Delete(&inIShortcut,rShortcut);
559 rDelete(rShortcut);
560 }
561 return std::pair<poly,int>(monomial,-1);
562}
int k
Definition: cfEzgcd.cc:99
g
Definition: cfModGcd.cc:4090
bool isValuationNonTrivial() const
ring getShortcutRingPrependingWeight(const ring r, const gfan::ZVector &w) const
If valuation trivial, returns a copy of r with a positive weight prepended, such that any ideal homog...
static FORCE_INLINE BOOLEAN n_IsOne(number n, const coeffs r)
TRUE iff 'n' represents the one element.
Definition: coeffs.h:465
poly searchForMonomialViaStepwiseSaturation(const ideal I, const ring r, const gfan::ZVector w0)
poly initial(const poly p, const ring r, const gfan::ZVector &w)
Returns the initial form of p with respect to w.
Definition: initial.cc:30
#define pNext(p)
Definition: monomials.h:36
#define p_GetCoeff(p, r)
Definition: monomials.h:50
poly p_PermPoly(poly p, const int *perm, const ring oldRing, const ring dst, nMapFunc nMap, const int *par_perm, int OldPar, BOOLEAN use_mult)
Definition: p_polys.cc:4126
static long p_GetExp(const poly p, const unsigned long iBitmask, const int VarOffset)
get a single variable exponent @Note: the integer VarOffset encodes:
Definition: p_polys.h:467
static short rVar(const ring r)
#define rVar(r) (r->N)
Definition: ring.h:592
ideal idInit(int idsize, int rank)
initialise an ideal / module
Definition: simpleideals.cc:35
void idSkipZeroes(ideal ide)
gives an ideal/module the minimal possible size
bool checkForNonPositiveEntries(const gfan::ZVector &w)

◆ computeFlip()

std::pair< ideal, ring > tropicalStrategy::computeFlip ( const ideal  Ir,
const ring  r,
const gfan::ZVector &  interiorPoint,
const gfan::ZVector &  facetNormal 
) const

given an interior point of a groebner cone computes the groebner cone adjacent to it

Definition at line 762 of file tropicalStrategy.cc.

765{
766 assume(isValuationTrivial() || interiorPoint[0].sign()<0);
768 assume(checkWeightVector(Ir,r,interiorPoint));
769
770 // get a generating system of the initial ideal
771 // and compute a standard basis with respect to adjacent ordering
772 ideal inIr = initial(Ir,r,interiorPoint);
773 ring sAdjusted = copyAndChangeOrderingWP(r,interiorPoint,facetNormal);
774 nMapFunc identity = n_SetMap(r->cf,sAdjusted->cf);
775 int k = IDELEMS(Ir);
776 ideal inIsAdjusted = idInit(k);
777 for (int i=0; i<k; i++)
778 inIsAdjusted->m[i] = p_PermPoly(inIr->m[i],NULL,r,sAdjusted,identity,NULL,0);
779 ideal inJsAdjusted = computeStdOfInitialIdeal(inIsAdjusted,sAdjusted);
780
781 // find witnesses of the new standard basis elements of the initial ideal
782 // with the help of the old standard basis of the ideal
783 k = IDELEMS(inJsAdjusted);
784 ideal inJr = idInit(k);
785 identity = n_SetMap(sAdjusted->cf,r->cf);
786 for (int i=0; i<k; i++)
787 inJr->m[i] = p_PermPoly(inJsAdjusted->m[i],NULL,sAdjusted,r,identity,NULL,0);
788
789 ideal Jr = computeWitness(inJr,inIr,Ir,r);
790 ring s = copyAndChangeOrderingLS(r,interiorPoint,facetNormal);
791 identity = n_SetMap(r->cf,s->cf);
792 ideal Js = idInit(k);
793 for (int i=0; i<k; i++)
794 Js->m[i] = p_PermPoly(Jr->m[i],NULL,r,s,identity,NULL,0);
795
796 reduce(Js,s);
797 assume(areIdealsEqual(Js,s,Ir,r));
799 assume(checkWeightVector(Js,s,interiorPoint));
800
801 // cleanup
802 id_Delete(&inIsAdjusted,sAdjusted);
803 id_Delete(&inJsAdjusted,sAdjusted);
804 rDelete(sAdjusted);
805 id_Delete(&inIr,r);
806 id_Delete(&Jr,r);
807 id_Delete(&inJr,r);
808
810 return std::make_pair(Js,s);
811}
ring copyAndChangeOrderingLS(const ring r, const gfan::ZVector &w, const gfan::ZVector &v) const
ideal computeWitness(const ideal inJ, const ideal inI, const ideal I, const ring r) const
suppose w a weight in maximal groebner cone of > suppose I (initially) reduced standard basis w....
ring copyAndChangeOrderingWP(const ring r, const gfan::ZVector &w, const gfan::ZVector &v) const
ideal computeStdOfInitialIdeal(const ideal inI, const ring r) const
given generators of the initial ideal, computes its standard basis
bool checkForUniformizingBinomial(const ideal I, const ring r) const
if valuation non-trivial, checks whether the generating system contains p-t otherwise returns true
bool isOrderingLocalInT(const ring r)
static int sign(int x)
Definition: ring.cc:3427
bool checkWeightVector(const ideal I, const ring r, const gfan::ZVector &weightVector, bool checkBorder)
bool areIdealsEqual(ideal I, ring r, ideal J, ring s)

◆ computeLift()

ideal tropicalStrategy::computeLift ( const ideal  inJs,
const ring  s,
const ideal  inIr,
const ideal  Ir,
const ring  r 
) const

Definition at line 688 of file tropicalStrategy.cc.

689{
690 int k = IDELEMS(inJs);
691 ideal inJr = idInit(k);
692 nMapFunc identitysr = n_SetMap(s->cf,r->cf);
693 for (int i=0; i<k; i++)
694 inJr->m[i] = p_PermPoly(inJs->m[i],NULL,s,r,identitysr,NULL,0);
695
696 ideal Jr = computeWitness(inJr,inIr,Ir,r);
697 nMapFunc identityrs = n_SetMap(r->cf,s->cf);
698 ideal Js = idInit(k);
699 for (int i=0; i<k; i++)
700 Js->m[i] = p_PermPoly(Jr->m[i],NULL,r,s,identityrs,NULL,0);
701 return Js;
702}

◆ computeStdOfInitialIdeal()

ideal tropicalStrategy::computeStdOfInitialIdeal ( const ideal  inI,
const ring  r 
) const

given generators of the initial ideal, computes its standard basis

Definition at line 656 of file tropicalStrategy.cc.

657{
658 // if valuation trivial, then compute std as usual
659 if (isValuationTrivial())
660 return gfanlib_kStd_wrapper(inI,r);
661
662 // if valuation non-trivial, then uniformizing parameter is in ideal
663 // so switch to residue field first and compute standard basis over the residue field
664 ring rShortcut = copyAndChangeCoefficientRing(r);
665 nMapFunc takingResidues = n_SetMap(r->cf,rShortcut->cf);
666 int k = IDELEMS(inI);
667 ideal inIShortcut = idInit(k);
668 for (int i=0; i<k; i++)
669 inIShortcut->m[i] = p_PermPoly(inI->m[i],NULL,r,rShortcut,takingResidues,NULL,0);
670 ideal inJShortcut = gfanlib_kStd_wrapper(inIShortcut,rShortcut);
671
672 // and lift the result back to the ring with valuation
673 nMapFunc takingRepresentatives = n_SetMap(rShortcut->cf,r->cf);
674 k = IDELEMS(inJShortcut);
675 ideal inJ = idInit(k+1);
676 inJ->m[0] = p_One(r);
677 nMapFunc identity = n_SetMap(startingRing->cf,r->cf);
678 p_SetCoeff(inJ->m[0],identity(uniformizingParameter,startingRing->cf,r->cf),r);
679 for (int i=0; i<k; i++)
680 inJ->m[i+1] = p_PermPoly(inJShortcut->m[i],NULL,rShortcut,r,takingRepresentatives,NULL,0);
681
682 id_Delete(&inJShortcut,rShortcut);
683 id_Delete(&inIShortcut,rShortcut);
684 rDelete(rShortcut);
685 return inJ;
686}
ring copyAndChangeCoefficientRing(const ring r) const
ideal gfanlib_kStd_wrapper(ideal I, ring r, tHomog h=testHomog)
Definition: std_wrapper.cc:6

◆ computeWitness()

ideal tropicalStrategy::computeWitness ( const ideal  inJ,
const ideal  inI,
const ideal  I,
const ring  r 
) const

suppose w a weight in maximal groebner cone of > suppose I (initially) reduced standard basis w.r.t.

> and inI initial forms of its elements w.r.t. w suppose inJ elements of initial ideal that are homogeneous w.r.t w returns J elements of ideal whose initial form w.r.t. w are inI in particular, if w lies also inthe maximal groebner cone of another ordering >' and inJ is a standard basis of the initial ideal w.r.t. >' then the returned J will be a standard baiss of the ideal w.r.t. >'

change ground ring into finite field and map the data into it

Compute a division with remainder over the finite field and map the result back to r

Compute the normal forms

Definition at line 574 of file tropicalStrategy.cc.

575{
576 // if the valuation is trivial and the ring and ideal have not been extended,
577 // then it is sufficient to return the difference between the elements of inJ
578 // and their normal forms with respect to I and r
579 if (isValuationTrivial())
580 return witness(inJ,I,r);
581 // if the valuation is non-trivial and the ring and ideal have been extended,
582 // then we can make a shortcut through the residue field
583 else
584 {
585 assume(IDELEMS(inI)==IDELEMS(I));
587 assume(uni>=0);
588 /**
589 * change ground ring into finite field
590 * and map the data into it
591 */
592 ring rShortcut = copyAndChangeCoefficientRing(r);
593
594 int k = IDELEMS(inJ);
595 int l = IDELEMS(I);
596 ideal inJShortcut = idInit(k);
597 ideal inIShortcut = idInit(l);
598 nMapFunc takingResidues = n_SetMap(r->cf,rShortcut->cf);
599 for (int i=0; i<k; i++)
600 inJShortcut->m[i] = p_PermPoly(inJ->m[i],NULL,r,rShortcut,takingResidues,NULL,0);
601 for (int j=0; j<l; j++)
602 inIShortcut->m[j] = p_PermPoly(inI->m[j],NULL,r,rShortcut,takingResidues,NULL,0);
603 id_Test(inJShortcut,rShortcut);
604 id_Test(inIShortcut,rShortcut);
605
606 /**
607 * Compute a division with remainder over the finite field
608 * and map the result back to r
609 */
610 matrix QShortcut = divisionDiscardingRemainder(inJShortcut,inIShortcut,rShortcut);
611 matrix Q = mpNew(l,k);
612 nMapFunc takingRepresentatives = n_SetMap(rShortcut->cf,r->cf);
613 for (int ij=k*l-1; ij>=0; ij--)
614 Q->m[ij] = p_PermPoly(QShortcut->m[ij],NULL,rShortcut,r,takingRepresentatives,NULL,0);
615
616 nMapFunc identity = n_SetMap(startingRing->cf,r->cf);
617 number p = identity(uniformizingParameter,startingRing->cf,r->cf);
618
619 /**
620 * Compute the normal forms
621 */
622 ideal J = idInit(k);
623 for (int j=0; j<k; j++)
624 {
625 poly q0 = p_Copy(inJ->m[j],r);
626 for (int i=0; i<l; i++)
627 {
628 poly qij = p_Copy(MATELEM(Q,i+1,j+1),r);
629 poly inIi = p_Copy(inI->m[i],r);
630 q0 = p_Add_q(q0,p_Neg(p_Mult_q(qij,inIi,r),r),r);
631 }
632 q0 = p_Div_nn(q0,p,r);
633 poly q0g0 = p_Mult_q(q0,p_Copy(I->m[uni],r),r);
634 // q0 = NULL;
635 poly qigi = NULL;
636 for (int i=0; i<l; i++)
637 {
638 poly qij = p_Copy(MATELEM(Q,i+1,j+1),r);
639 // poly inIi = p_Copy(I->m[i],r);
640 poly Ii = p_Copy(I->m[i],r);
641 qigi = p_Add_q(qigi,p_Mult_q(qij,Ii,r),r);
642 }
643 J->m[j] = p_Add_q(q0g0,qigi,r);
644 }
645
646 id_Delete(&inIShortcut,rShortcut);
647 id_Delete(&inJShortcut,rShortcut);
648 mp_Delete(&QShortcut,rShortcut);
649 rDelete(rShortcut);
650 mp_Delete(&Q,r);
651 n_Delete(&p,r->cf);
652 return J;
653 }
654}
int l
Definition: cfEzgcd.cc:100
poly * m
Definition: matpol.h:18
int findPositionOfUniformizingBinomial(const ideal I, const ring r) const
int j
Definition: facHensel.cc:110
STATIC_VAR jList * Q
Definition: janet.cc:30
void mp_Delete(matrix *a, const ring r)
Definition: matpol.cc:880
matrix mpNew(int r, int c)
create a r x c zero-matrix
Definition: matpol.cc:37
#define MATELEM(mat, i, j)
1-based access to matrix
Definition: matpol.h:29
poly p_Div_nn(poly p, const number n, const ring r)
Definition: p_polys.cc:1501
static poly p_Mult_q(poly p, poly q, const ring r)
Definition: p_polys.h:1112
static poly p_Copy(poly p, const ring r)
returns a copy of p
Definition: p_polys.h:844
matrix divisionDiscardingRemainder(const poly f, const ideal G, const ring r)
Computes a division discarding remainder of f with respect to G.
Definition: witness.cc:9
poly witness(const poly m, const ideal I, const ideal inI, const ring r)
Let w be the uppermost weight vector in the matrix defining the ordering on r.
Definition: witness.cc:34

◆ copyAndChangeCoefficientRing()

ring tropicalStrategy::copyAndChangeCoefficientRing ( const ring  r) const
private

Definition at line 564 of file tropicalStrategy.cc.

565{
566 ring rShortcut = rCopy0(r,FALSE); // do not copy q-ideal
567 nKillChar(rShortcut->cf);
568 rShortcut->cf = nCopyCoeff(shortcutRing->cf);
569 rComplete(rShortcut);
570 rTest(rShortcut);
571 return rShortcut;
572}
static FORCE_INLINE coeffs nCopyCoeff(const coeffs r)
"copy" coeffs, i.e. increment ref
Definition: coeffs.h:430

◆ copyAndChangeOrderingLS()

ring tropicalStrategy::copyAndChangeOrderingLS ( const ring  r,
const gfan::ZVector &  w,
const gfan::ZVector &  v 
) const
private

Definition at line 734 of file tropicalStrategy.cc.

735{
736 // copy shortcutRing and change to desired ordering
737 bool ok;
738 ring s = rCopy0(r,FALSE,FALSE);
739 int n = rVar(s);
740 s->order = (rRingOrder_t*) omAlloc0(5*sizeof(rRingOrder_t));
741 s->block0 = (int*) omAlloc0(5*sizeof(int));
742 s->block1 = (int*) omAlloc0(5*sizeof(int));
743 s->wvhdl = (int**) omAlloc0(5*sizeof(int**));
744 s->order[0] = ringorder_a;
745 s->block0[0] = 1;
746 s->block1[0] = n;
747 s->wvhdl[0] = ZVectorToIntStar(w,ok);
748 s->order[1] = ringorder_a;
749 s->block0[1] = 1;
750 s->block1[1] = n;
751 s->wvhdl[1] = ZVectorToIntStar(v,ok);
752 s->order[2] = ringorder_lp;
753 s->block0[2] = 1;
754 s->block1[2] = n;
755 s->order[3] = ringorder_C;
756 rComplete(s);
757 rTest(s);
758
759 return s;
760}
int * ZVectorToIntStar(const gfan::ZVector &v, bool &overflow)
#define omAlloc0(size)
Definition: omAllocDecl.h:211
rRingOrder_t
order stuff
Definition: ring.h:68
@ ringorder_lp
Definition: ring.h:77
@ ringorder_a
Definition: ring.h:70
@ ringorder_C
Definition: ring.h:73

◆ copyAndChangeOrderingWP()

ring tropicalStrategy::copyAndChangeOrderingWP ( const ring  r,
const gfan::ZVector &  w,
const gfan::ZVector &  v 
) const
private

Definition at line 704 of file tropicalStrategy.cc.

705{
706 // copy shortcutRing and change to desired ordering
707 bool ok;
708 ring s = rCopy0(r,FALSE,FALSE);
709 int n = rVar(s);
710 gfan::ZVector wAdjusted = adjustWeightForHomogeneity(w);
711 gfan::ZVector vAdjusted = adjustWeightUnderHomogeneity(v,wAdjusted);
712 s->order = (rRingOrder_t*) omAlloc0(5*sizeof(rRingOrder_t));
713 s->block0 = (int*) omAlloc0(5*sizeof(int));
714 s->block1 = (int*) omAlloc0(5*sizeof(int));
715 s->wvhdl = (int**) omAlloc0(5*sizeof(int**));
716 s->order[0] = ringorder_a;
717 s->block0[0] = 1;
718 s->block1[0] = n;
719 s->wvhdl[0] = ZVectorToIntStar(wAdjusted,ok);
720 s->order[1] = ringorder_a;
721 s->block0[1] = 1;
722 s->block1[1] = n;
723 s->wvhdl[1] = ZVectorToIntStar(vAdjusted,ok);
724 s->order[2] = ringorder_lp;
725 s->block0[2] = 1;
726 s->block1[2] = n;
727 s->order[3] = ringorder_C;
728 rComplete(s);
729 rTest(s);
730
731 return s;
732}
gfan::ZVector adjustWeightUnderHomogeneity(gfan::ZVector v, gfan::ZVector w) const
Given strictly positive weight w and weight v, returns a strictly positive weight u such that on an i...
gfan::ZVector adjustWeightForHomogeneity(gfan::ZVector w) const
Given weight w, returns a strictly positive weight u such that an ideal satisfying the valuation-sepc...

◆ findPositionOfUniformizingBinomial()

int tropicalStrategy::findPositionOfUniformizingBinomial ( const ideal  I,
const ring  r 
) const
private

Definition at line 842 of file tropicalStrategy.cc.

843{
845
846 // if the valuation is non-trivial then checks if the first generator is p-t
847 nMapFunc identity = n_SetMap(startingRing->cf,r->cf);
848 poly p = p_One(r);
849 p_SetCoeff(p,identity(uniformizingParameter,startingRing->cf,r->cf),r);
850 poly t = p_One(r);
851 p_SetExp(t,1,1,r);
852 p_Setm(t,r);
853 poly pt = p_Add_q(p,p_Neg(t,r),r);
854
855 for (int i=0; i<IDELEMS(I); i++)
856 {
857 if (p_EqualPolys(I->m[i],pt,r))
858 {
859 p_Delete(&pt,r);
860 return i;
861 }
862 }
863 p_Delete(&pt,r);
864 return -1;
865}

◆ getExpectedAmbientDimension()

int tropicalStrategy::getExpectedAmbientDimension ( ) const
inline

Definition at line 191 of file tropicalStrategy.h.

192 {
193 return rVar(startingRing);
194 }

◆ getExpectedDimension()

int tropicalStrategy::getExpectedDimension ( ) const
inline

returns the expected Dimension of the polyhedral output

Definition at line 199 of file tropicalStrategy.h.

200 {
201 return expectedDimension;
202 }

◆ getHomogeneitySpace()

gfan::ZCone tropicalStrategy::getHomogeneitySpace ( ) const
inline

returns the homogeneity space of the preimage ideal

Definition at line 222 of file tropicalStrategy.h.

223 {
224 return linealitySpace;
225 }

◆ getOriginalIdeal()

ideal tropicalStrategy::getOriginalIdeal ( ) const
inline

returns the input ideal over the field with valuation

Definition at line 167 of file tropicalStrategy.h.

168 {
170 return originalIdeal;
171 }

◆ getOriginalRing()

ring tropicalStrategy::getOriginalRing ( ) const
inline

returns the polynomial ring over the field with valuation

Definition at line 158 of file tropicalStrategy.h.

159 {
161 return originalRing;
162 }

◆ getShortcutRing()

ring tropicalStrategy::getShortcutRing ( ) const
inline

Definition at line 213 of file tropicalStrategy.h.

214 {
216 return shortcutRing;
217 }

◆ getShortcutRingPrependingWeight()

ring tropicalStrategy::getShortcutRingPrependingWeight ( const ring  r,
const gfan::ZVector &  w 
) const

If valuation trivial, returns a copy of r with a positive weight prepended, such that any ideal homogeneous with respect to w is homogeneous with respect to that weight.

If valuation non-trivial, changes the coefficient ring to the residue field.

Definition at line 448 of file tropicalStrategy.cc.

449{
450 ring rShortcut = rCopy0(r,FALSE); // do not copy q-ideal
451
452 // save old ordering
453 rRingOrder_t* order = rShortcut->order;
454 int* block0 = rShortcut->block0;
455 int* block1 = rShortcut->block1;
456 int** wvhdl = rShortcut->wvhdl;
457
458 // adjust weight and create new ordering
459 gfan::ZVector w = adjustWeightForHomogeneity(v);
460 int h = rBlocks(r); int n = rVar(r);
461 rShortcut->order = (rRingOrder_t*) omAlloc0((h+2)*sizeof(rRingOrder_t));
462 rShortcut->block0 = (int*) omAlloc0((h+2)*sizeof(int));
463 rShortcut->block1 = (int*) omAlloc0((h+2)*sizeof(int));
464 rShortcut->wvhdl = (int**) omAlloc0((h+2)*sizeof(int*));
465 rShortcut->order[0] = ringorder_a;
466 rShortcut->block0[0] = 1;
467 rShortcut->block1[0] = n;
468 bool overflow;
469 rShortcut->wvhdl[0] = ZVectorToIntStar(w,overflow);
470 for (int i=1; i<=h; i++)
471 {
472 rShortcut->order[i] = order[i-1];
473 rShortcut->block0[i] = block0[i-1];
474 rShortcut->block1[i] = block1[i-1];
475 rShortcut->wvhdl[i] = wvhdl[i-1];
476 }
477
478 // if valuation non-trivial, change coefficient ring to residue field
480 {
481 nKillChar(rShortcut->cf);
482 rShortcut->cf = nCopyCoeff(shortcutRing->cf);
483 }
484 rComplete(rShortcut);
485 rTest(rShortcut);
486
487 // delete old ordering
488 omFree(order);
489 omFree(block0);
490 omFree(block1);
491 omFree(wvhdl);
492
493 return rShortcut;
494}
STATIC_VAR Poly * h
Definition: janet.cc:971
#define omFree(addr)
Definition: omAllocDecl.h:261
static int rBlocks(const ring r)
Definition: ring.h:568

◆ getStartingIdeal()

ideal tropicalStrategy::getStartingIdeal ( ) const
inline

returns the input ideal

Definition at line 185 of file tropicalStrategy.h.

186 {
188 return startingIdeal;
189 }

◆ getStartingRing()

ring tropicalStrategy::getStartingRing ( ) const
inline

returns the polynomial ring over the valuation ring

Definition at line 176 of file tropicalStrategy.h.

177 {
179 return startingRing;
180 }

◆ getUniformizingParameter()

number tropicalStrategy::getUniformizingParameter ( ) const
inline

returns the uniformizing parameter in the valuation ring

Definition at line 207 of file tropicalStrategy.h.

◆ homogeneitySpaceContains()

bool tropicalStrategy::homogeneitySpaceContains ( const gfan::ZVector &  v) const
inline

returns true, if v is contained in the homogeneity space; false otherwise

Definition at line 230 of file tropicalStrategy.h.

231 {
232 return linealitySpace.contains(v);
233 }

◆ isValuationNonTrivial()

bool tropicalStrategy::isValuationNonTrivial ( ) const
inline

Definition at line 149 of file tropicalStrategy.h.

150 {
151 bool b = (uniformizingParameter!=NULL);
152 return b;
153 }
CanonicalForm b
Definition: cfModGcd.cc:4103

◆ isValuationTrivial()

bool tropicalStrategy::isValuationTrivial ( ) const
inline

Definition at line 144 of file tropicalStrategy.h.

145 {
146 bool b = (uniformizingParameter==NULL);
147 return b;
148 }

◆ negateWeight()

gfan::ZVector tropicalStrategy::negateWeight ( const gfan::ZVector &  w) const
inline

Definition at line 263 of file tropicalStrategy.h.

264 {
265 gfan::ZVector wNeg(w.size());
266
267 if (this->isValuationNonTrivial())
268 {
269 wNeg[0]=w[0];
270 for (unsigned i=1; i<w.size(); i++)
271 wNeg[i]=w[i];
272 }
273 else
274 wNeg = -w;
275
276 return wNeg;
277 }

◆ operator=()

tropicalStrategy & tropicalStrategy::operator= ( const tropicalStrategy currentStrategy)

assignment operator

Definition at line 363 of file tropicalStrategy.cc.

◆ pReduce()

void tropicalStrategy::pReduce ( ideal  I,
const ring  r 
) const

Definition at line 432 of file tropicalStrategy.cc.

433{
434 rTest(r);
435 id_Test(I,r);
436
437 if (isValuationTrivial())
438 return;
439
440 nMapFunc identity = n_SetMap(startingRing->cf,r->cf);
441 number p = identity(uniformizingParameter,startingRing->cf,r->cf);
442 ::pReduce(I,p,r);
443 n_Delete(&p,r->cf);
444
445 return;
446}
void pReduce(ideal I, const ring r) const

◆ putUniformizingBinomialInFront()

void tropicalStrategy::putUniformizingBinomialInFront ( ideal  I,
const ring  r,
const number  q 
) const
private

Definition at line 387 of file tropicalStrategy.cc.

388{
389 poly p = p_One(r);
390 p_SetCoeff(p,q,r);
391 poly t = p_One(r);
392 p_SetExp(t,1,1,r);
393 p_Setm(t,r);
394 poly pt = p_Add_q(p,p_Neg(t,r),r);
395
396 int k = IDELEMS(I);
397 int l;
398 for (l=0; l<k; l++)
399 {
400 if (p_EqualPolys(I->m[l],pt,r))
401 break;
402 }
403 p_Delete(&pt,r);
404
405 if (l>1)
406 {
407 pt = I->m[l];
408 for (int i=l; i>0; i--)
409 I->m[l] = I->m[l-1];
410 I->m[0] = pt;
411 pt = NULL;
412 }
413 return;
414}

◆ reduce()

bool tropicalStrategy::reduce ( ideal  I,
const ring  r 
) const

reduces the generators of an ideal I so that the inequalities and equations of the Groebner cone can be read off.

Definition at line 416 of file tropicalStrategy.cc.

417{
418 rTest(r);
419 id_Test(I,r);
420
421 nMapFunc identity = n_SetMap(startingRing->cf,r->cf);
422 number p = NULL;
424 p = identity(uniformizingParameter,startingRing->cf,r->cf);
425 bool b = extraReductionAlgorithm(I,r,p);
426 // putUniformizingBinomialInFront(I,r,p);
427 if (p!=NULL) n_Delete(&p,r->cf);
428
429 return b;
430}

◆ restrictToLowerHalfSpace()

bool tropicalStrategy::restrictToLowerHalfSpace ( ) const
inline

returns true, if valuation non-trivial, false otherwise

Definition at line 238 of file tropicalStrategy.h.

239 {
240 return onlyLowerHalfSpace;
241 }

Field Documentation

◆ expectedDimension

int tropicalStrategy::expectedDimension
private

the expected Dimension of the polyhedral output, i.e.

the dimension of the ideal if valuation trivial or the dimension of the ideal plus one if valuation non-trivial (as the output is supposed to be intersected with a hyperplane)

Definition at line 53 of file tropicalStrategy.h.

◆ extraReductionAlgorithm

bool(* tropicalStrategy::extraReductionAlgorithm) (ideal I, ring r, number p)
private

A function that reduces the generators of an ideal I so that the inequalities and equations of the Groebner cone can be read off.

Definition at line 98 of file tropicalStrategy.h.

◆ linealitySpace

gfan::ZCone tropicalStrategy::linealitySpace
private

the homogeneity space of the Grobner fan

Definition at line 57 of file tropicalStrategy.h.

◆ onlyLowerHalfSpace

bool tropicalStrategy::onlyLowerHalfSpace
private

true if valuation non-trivial, false otherwise

Definition at line 78 of file tropicalStrategy.h.

◆ originalIdeal

ideal tropicalStrategy::originalIdeal
private

input ideal, assumed to be a homogeneous prime ideal

Definition at line 46 of file tropicalStrategy.h.

◆ originalRing

ring tropicalStrategy::originalRing
private

polynomial ring over a field with valuation

Definition at line 42 of file tropicalStrategy.h.

◆ shortcutRing

ring tropicalStrategy::shortcutRing
private

polynomial ring over the residue field

Definition at line 73 of file tropicalStrategy.h.

◆ startingIdeal

ideal tropicalStrategy::startingIdeal
private

preimage of the input ideal under the map that sends t to the uniformizing parameter

Definition at line 65 of file tropicalStrategy.h.

◆ startingRing

ring tropicalStrategy::startingRing
private

polynomial ring over the valuation ring extended by one extra variable t

Definition at line 61 of file tropicalStrategy.h.

◆ uniformizingParameter

number tropicalStrategy::uniformizingParameter
private

uniformizing parameter in the valuation ring

Definition at line 69 of file tropicalStrategy.h.

◆ weightAdjustingAlgorithm1

gfan::ZVector(* tropicalStrategy::weightAdjustingAlgorithm1) (const gfan::ZVector &w)
private

A function such that: Given weight w, returns a strictly positive weight u such that an ideal satisfying the valuation-sepcific homogeneity conditions is weighted homogeneous with respect to w if and only if it is homogeneous with respect to u.

Definition at line 86 of file tropicalStrategy.h.

◆ weightAdjustingAlgorithm2

gfan::ZVector(* tropicalStrategy::weightAdjustingAlgorithm2) (const gfan::ZVector &v, const gfan::ZVector &w)
private

A function such that: Given strictly positive weight w and weight v, returns a strictly positive weight u such that on an ideal that is weighted homogeneous with respect to w the weights u and v coincide.

Definition at line 93 of file tropicalStrategy.h.


The documentation for this class was generated from the following files: