FeynCalc / feyncalc

Mathematica package for algebraic calculations in elementary particle physics.
https://feyncalc.github.io
GNU General Public License v3.0
309 stars 87 forks source link

Loop integral with FCMultiLoopTID #273

Closed liandk closed 1 week ago

liandk commented 2 months ago
FAD[{q - p - k, 0, 1}, {q - p, 0, 1}, {p, 0, 1}, {k, 0, 1}] FVD[
    k, \[Mu]] TR[
    GSD[q - p - k] . GAD[\[Rho]] . GSD[q - p] . GAD[\[Alpha]] . 
     GSD[p] . GAD[\[Mu]]] // Calc;
FCMultiLoopTID[%, {p, k}] // Calc;
ToFI[%, {k, p}, {q}] // Calc;
TarcerRecurse[%]

The result is 屏幕截图 2024-06-24 182828

liandk commented 1 month ago

Sorry, this issue is duplicate. So I close it.

vsht commented 1 month ago

It's because if you do the shift, then after FCMultiLoopTID one of the integrals (FAD[q, q - p]) doesn't depend on k but is not automatically set to zero.

However, filtering out such integrals in full generality is quite tricky because also cases such as FAD[q+k, q + k - p] could appear: it looks like a perfectly fine 2L integral, but is actually only 1L.

For now I would recommend you to use the approach from the FeynCalc 10 examples with FCLoopFindTopologies and so on.

FCMultiLoopTID is not a very clever or robust function and was only meant as a placeholder in version 9 until something better comes up. Perhaps it would even make sense to get rid of it completely in the future versions.

liandk commented 1 month ago

Thanks for your answer.

vsht commented 1 month ago

Here's a sample code for that. The new approach is much better than FCMultiLoopTID and all other old tricks to deal with multiloop integrals in FeynCalc

$LoadAddOns = {"FeynHelpers"};

<< FeynCalc`

FCClearScalarProducts[];
ScalarProduct[q, q] = qq;

ampRaw1 = 
 FAD[{q - p - k, 0, 1}, {q - p, 0, 1}, {p, 0, 1}, {k, 0, 1}] FVD[
   k, \[Mu]] TR[
   GSD[q - p - k] . GAD[\[Rho]] . GSD[q - p] . GAD[\[Alpha]] . 
    GSD[p] . GAD[\[Mu]]]

{amp, topos} = 
 FCLoopFindTopologies[ampRaw1 /. {p -> q - k - p}, {p, k}]

toposC = FCLoopBasisFindCompletion[topos]
subtoposC = FCLoopFindSubtopologies[toposC];

mappings = 
 FCLoopFindTopologyMappings[topos, PreferredTopologies -> subtoposC]

ampReduced = FCLoopTensorReduce[amp, topos];

ampPreFinal = FCLoopApplyTopologyMappings[ampReduced, mappings];

ampPreFinal

FIREPrepareStartFile[mappings[[2]], FCGetNotebookDirectory[]]

FIRECreateStartFile[FCGetNotebookDirectory[], mappings[[2]]]

FIRECreateConfigFile[mappings[[2]], FCGetNotebookDirectory[]]

FIRECreateIntegralFile[Cases2[ampPreFinal, GLI], mappings[[2]], 
 FCGetNotebookDirectory[]]

FIRERunReduction[FCGetNotebookDirectory[], mappings[[2]]]

tables = 
  FIREImportResults[mappings[[2]], FCGetNotebookDirectory[]] // 
   Flatten;

resPreFinal = Collect2[ampPreFinal /. tables, GLI]

Bildschirmfoto vom 2024-07-01 11-31-39

vsht commented 2 weeks ago

@liandk I added a fix that should prevent such issues in the future, even though FDS might get a bit slower now. But it's a fair price to pay for a guarantee that all scaleless integrals are being set to zero early on.

Can you please check if the dev version now works for all your cases?