JHUGen MELA  v2.4.1
Matrix element calculations as used in JHUGen. MELA is an important tool that was used for the Higgs boson discovery and for precise measurements of its structure and interactions. Please see the website https://spin.pha.jhu.edu/ and papers cited there for more details, and kindly cite those papers when using this code.
Functions/Subroutines
buildtensors Module Reference

Functions/Subroutines

subroutine calctensora_list (TA, TAuv, TAerr, CoefsA, CoefsAuv, CoefsAerr, rmax)
 
subroutine calctensora (TA, TAuv, TAerr, CoefsA, CoefsAuv, CoefsAerr, rmax)
 
subroutine calctensorb_list (TB, TBuv, TBerr, CoefsB, CoefsBuv, CoefsBerr, mom, rmax)
 
subroutine calctensorb (TB, TBuv, TBerr, CoefsB, CoefsBuv, CoefsBerr, mom, rmax)
 
subroutine calctensorc_list (TC, TCuv, TCerr, CoefsC, CoefsCuv, CoefsCerr, MomVec, rmax)
 
subroutine calctensorc (TC, TCuv, TCerr, CoefsC, CoefsCuv, CoefsCerr, MomVec, rmax)
 
subroutine calctensord_list (TD, TDuv, TDerr, CoefsD, CoefsDuv, CoefsDerr, MomVec, rmax)
 
subroutine calctensord (TD, TDuv, TDerr, CoefsD, CoefsDuv, CoefsDerr, MomVec, rmax)
 
subroutine calctensore_list (TE, TEuv, TEerr, CoefsE, CoefsEuv, CoefsEerr, MomVec, rmax)
 
subroutine calctensore (TE, TEuv, TEerr, CoefsE, CoefsEuv, CoefsEerr, MomVec, rmax)
 
subroutine calctensorf_list (TF, TFuv, TFerr, CoefsF, CoefsFuv, CoefsFerr, MomVec, rmax)
 
subroutine calctensorf (TF, TFuv, TFerr, CoefsF, CoefsFuv, CoefsFerr, MomVec, rmax)
 
subroutine calctensorfuv_list (TFuv, CoefsFuv, MomVec, rmax)
 
subroutine calctensorfuv (TFuv, CoefsFuv, MomVec, rmax)
 
subroutine calctensorg_list (TG, TGuv, TGerr, CoefsG, CoefsGuv, CoefsGerr, MomVec, rmax)
 
subroutine calctensorg (TG, TGuv, TGerr, CoefsG, CoefsGuv, CoefsGerr, MomVec, rmax)
 
subroutine calctensortn_list (TN, TNuv, TNerr, CoefsN, CoefsNuv, CoefsNerr, MomVec, N, rmax)
 
subroutine calctensortn (TN, TNuv, TNerr, CoefsN, CoefsNuv, CoefsNerr, MomVec, N, rmax)
 
subroutine calctensortnuv_list (TNuv, CoefsNuv, MomVec, N, rmax)
 

Function/Subroutine Documentation

◆ calctensora()

subroutine buildtensors::calctensora ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TA,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TAuv,
double precision, dimension(0:rmax), intent(out)  TAerr,
double complex, dimension(0:rmax/2), intent(in)  CoefsA,
double complex, dimension(0:rmax/2), intent(in)  CoefsAuv,
double precision, dimension(0:rmax), intent(in)  CoefsAerr,
integer, intent(in)  rmax 
)

Definition at line 89 of file BuildTensors.F90.

89 
90  integer, intent(in) :: rmax
91  double complex, intent(in) :: CoefsA(0:rmax/2), CoefsAuv(0:rmax/2)
92  double precision, intent(in) :: CoefsAerr(0:rmax)
93  double complex, intent(out) :: TA(0:rmax,0:rmax,0:rmax,0:rmax)
94  double complex, intent(out) :: TAuv(0:rmax,0:rmax,0:rmax,0:rmax)
95  double precision, intent(out) :: TAerr(0:rmax)
96  double complex :: TA_aux(RtS(rmax)), TAuv_aux(RtS(rmax))
97  integer :: mu
98 
99  call calctensora_list(ta_aux,tauv_aux,taerr,coefsa,coefsauv,coefsaerr,rmax)
100 
101  do mu=1,rts(rmax)
102  ta(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = ta_aux(mu)
103  end do
104 
105  if (calcuv_ten) then
106  do mu=1,rts(rmax)
107  tauv(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = tauv_aux(mu)
108  end do
109  end if
110 

◆ calctensora_list()

subroutine buildtensors::calctensora_list ( double complex, dimension(rts(rmax)), intent(out)  TA,
double complex, dimension(rts(rmax)), intent(out)  TAuv,
double precision, dimension(0:rmax)  TAerr,
double complex, dimension(0:rmax/2), intent(in)  CoefsA,
double complex, dimension(0:rmax/2), intent(in)  CoefsAuv,
double precision, dimension(0:rmax), intent(in)  CoefsAerr,
integer, intent(in)  rmax 
)

Definition at line 43 of file BuildTensors.F90.

43 
44  integer, intent(in) :: rmax
45  double complex, intent(in) :: CoefsA(0:rmax/2), CoefsAuv(0:rmax/2)
46  double precision, intent(in) :: CoefsAerr(0:rmax)
47  double complex, intent(out) :: TA(RtS(rmax)), TAuv(RtS(rmax))
48  double precision :: TAerr(0:rmax)
49  double complex :: CA
50  integer :: mu1,mu2,nsum,mu,nu,r,a,cnt
51 
52  ta = 0d0
53  ta(1) = coefsa(0)
54  taerr = coefsaerr(0)
55 
56  do nsum=1,rmax/2
57  ca = coefsa(nsum)
58  do nu=rts(nsum-1)+1,rts(nsum)
59  ta(addgtab(1,nu)) = ta(addgtab(1,nu)) + ca*cftab(1,nu)
60  end do
61  end do
62 
63  if (calcuv_ten) then
64  tauv = 0d0
65  tauv(1) = coefsauv(0)
66 
67  do nsum=1,rmax/2
68  ca = coefsauv(nsum)
69  do nu=rts(nsum-1)+1,rts(nsum)
70  tauv(addgtab(1,nu)) = tauv(addgtab(1,nu)) + ca*cftab(1,nu)
71  end do
72  end do
73 
74  end if
75 
76 

◆ calctensorb()

subroutine buildtensors::calctensorb ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TB,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TBuv,
double precision, dimension(0:rmax), intent(out)  TBerr,
double complex, dimension(0:rmax/2,0:rmax), intent(in)  CoefsB,
double complex, dimension(0:rmax/2,0:rmax), intent(in)  CoefsBuv,
double precision, dimension(0:rmax), intent(in)  CoefsBerr,
double complex, dimension(0:3), intent(in)  mom,
integer, intent(in)  rmax 
)

Definition at line 229 of file BuildTensors.F90.

229 
230  integer, intent(in) :: rmax
231  double complex, intent(in) :: mom(0:3)
232  double complex, intent(in) :: CoefsB(0:rmax/2,0:rmax), CoefsBuv(0:rmax/2,0:rmax)
233  double precision, intent(in) :: CoefsBerr(0:rmax)
234  double complex, intent(out) :: TB(0:rmax,0:rmax,0:rmax,0:rmax)
235  double complex, intent(out) :: TBuv(0:rmax,0:rmax,0:rmax,0:rmax)
236  double precision, intent(out) :: TBerr(0:rmax)
237  double complex :: TB_aux(RtS(rmax)), TBuv_aux(RtS(rmax))
238  integer :: mu
239 
240  call calctensorb_list(tb_aux,tbuv_aux,tberr,coefsb,coefsbuv,coefsberr,mom,rmax)
241 
242  do mu=1,rts(rmax)
243  tb(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = tb_aux(mu)
244  end do
245 
246  if (calcuv_ten) then
247  do mu=1,rts(rmax)
248  tbuv(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = tbuv_aux(mu)
249  end do
250  end if
251 
252 

◆ calctensorb_list()

subroutine buildtensors::calctensorb_list ( double complex, dimension(rts(rmax)), intent(out)  TB,
double complex, dimension(rts(rmax)), intent(out)  TBuv,
double precision, dimension(0:rmax), intent(out)  TBerr,
double complex, dimension(0:rmax/2,0:rmax), intent(in)  CoefsB,
double complex, dimension(0:rmax/2,0:rmax), intent(in)  CoefsBuv,
double precision, dimension(0:rmax), intent(in)  CoefsBerr,
double complex, dimension(0:3), intent(in)  mom,
integer, intent(in)  rmax 
)

Definition at line 123 of file BuildTensors.F90.

123 
124  integer, intent(in) :: rmax
125  double complex, intent(in) :: mom(0:3)
126  double complex, intent(in) :: CoefsB(0:rmax/2,0:rmax), CoefsBuv(0:rmax/2,0:rmax)
127  double precision, intent(in) :: CoefsBerr(0:rmax)
128  double complex, intent(out) :: TB(RtS(rmax)), TBuv(RtS(rmax))
129  double precision, intent(out) :: TBerr(0:rmax)
130  double complex :: MomTen(RtS(rmax)), CB, Pmu
131  double precision :: MomMax
132  integer :: mu1,mu2,nsum,mu,nu,r,a,cnt
133 
134 
135  tb = 0d0
136  tb(1) = coefsb(0,0)
137  tberr = 0d0
138  tberr(0) = coefsberr(0)
139 
140  do nsum=1,rmax/2
141  cb = coefsb(nsum,0)
142 
143  do nu=rts(nsum-1)+1,rts(nsum)
144  tb(addgtab(1,nu)) = tb(addgtab(1,nu)) + cb*cftab(1,nu)
145  end do
146 
147  end do
148 
149  if (calcuv_ten) then
150  tbuv = 0d0
151  tbuv(1) = coefsbuv(0,0)
152 
153  do nsum=1,rmax/2
154  cb = coefsbuv(nsum,0)
155 
156  do nu=rts(nsum-1)+1,rts(nsum)
157  tbuv(addgtab(1,nu)) = tbuv(addgtab(1,nu)) + cb*cftab(1,nu)
158  end do
159 
160  end do
161  end if
162 
163  momten(1) = 1
164  do r=1,rmax
165 
166  mu1 = rts(r-1)+1
167  mu2 = rts(r)
168 
169  cnt = mu1
170  do mu=0,3
171  pmu = mom(mu)
172  do a = mu1-binomtable(r-1,r+2-mu),mu1-1
173  momten(cnt)=momten(a)*pmu
174  cnt = cnt+1
175  end do
176  end do
177 
178  cb = coefsb(0,r)
179  tb(mu1:mu2) = tb(mu1:mu2) + cb*momten(mu1:mu2)
180 
181  do nsum=1,(rmax-r)/2
182  cb = coefsb(nsum,r)
183 
184  do nu=rts(nsum-1)+1,rts(nsum)
185  do mu=mu1,mu2
186  tb(addgtab(mu,nu)) = tb(addgtab(mu,nu)) + cb*momten(mu)*cftab(mu,nu)
187  end do
188  end do
189 
190  end do
191 
192  if (calcuv_ten) then
193  cb = coefsbuv(0,r)
194  tbuv(mu1:mu2) = tbuv(mu1:mu2) + cb*momten(mu1:mu2)
195 
196  do nsum=1,(rmax-r)/2
197  cb = coefsbuv(nsum,r)
198 
199  do nu=rts(nsum-1)+1,rts(nsum)
200  do mu=mu1,mu2
201  tbuv(addgtab(mu,nu)) = tbuv(addgtab(mu,nu)) + cb*momten(mu)*cftab(mu,nu)
202  end do
203  end do
204 
205  end do
206 
207  end if
208 
209  end do
210 
211  mommax = maxval(abs(mom))
212  do r=1,rmax
213  tberr(r) = coefsberr(r)*mommax**r
214  end do
215 
216 

◆ calctensorc()

subroutine buildtensors::calctensorc ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TC,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TCuv,
double precision, dimension(0:rmax), intent(out)  TCerr,
double complex, dimension(0:rmax/2,0:rmax,0:rmax), intent(in)  CoefsC,
double complex, dimension(0:rmax/2,0:rmax,0:rmax), intent(in)  CoefsCuv,
double precision, dimension(0:rmax), intent(in)  CoefsCerr,
double complex, dimension(0:3,2), intent(in)  MomVec,
integer, intent(in)  rmax 
)

Definition at line 413 of file BuildTensors.F90.

413 
414  integer, intent(in) :: rmax
415  double complex, intent(in) :: MomVec(0:3,2)
416  double complex, intent(in) :: CoefsC(0:rmax/2,0:rmax,0:rmax)
417  double complex, intent(in) :: CoefsCuv(0:rmax/2,0:rmax,0:rmax)
418  double precision, intent(in) :: CoefsCerr(0:rmax)
419  double complex, intent(out) :: TC(0:rmax,0:rmax,0:rmax,0:rmax)
420  double complex, intent(out) :: TCuv(0:rmax,0:rmax,0:rmax,0:rmax)
421  double precision, intent(out) :: TCerr(0:rmax)
422  double complex :: TC_aux(RtS(rmax)), TCuv_aux(RtS(rmax))
423  integer :: mu
424 
425  call calctensorc_list(tc_aux,tcuv_aux,tcerr,coefsc,coefscuv,coefscerr,momvec,rmax)
426 
427  do mu=1,rts(rmax)
428  tc(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = tc_aux(mu)
429  end do
430 
431  if (calcuv_ten) then
432  do mu=1,rts(rmax)
433  tcuv(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = tcuv_aux(mu)
434  end do
435  else
436  tcuv = 0d0
437  end if
438 
439 

◆ calctensorc_list()

subroutine buildtensors::calctensorc_list ( double complex, dimension(rts(rmax)), intent(out)  TC,
double complex, dimension(rts(rmax)), intent(out)  TCuv,
double precision, dimension(0:rmax), intent(out)  TCerr,
double complex, dimension(0:rmax/2,0:rmax,0:rmax), intent(in)  CoefsC,
double complex, dimension(0:rmax/2,0:rmax,0:rmax), intent(in)  CoefsCuv,
double precision, dimension(0:rmax), intent(in)  CoefsCerr,
double complex, dimension(0:3,2), intent(in)  MomVec,
integer, intent(in)  rmax 
)

Definition at line 265 of file BuildTensors.F90.

265 
266  integer, intent(in) :: rmax
267  double complex, intent(in) :: MomVec(0:3,2)
268  double complex, intent(in) :: CoefsC(0:rmax/2,0:rmax,0:rmax)
269  double complex, intent(in) :: CoefsCuv(0:rmax/2,0:rmax,0:rmax)
270  double precision, intent(in) :: CoefsCerr(0:rmax)
271  double complex, intent(out) :: TC(RtS(rmax)), TCuv(RtS(rmax))
272  double precision, intent(out) :: TCerr(0:rmax)
273  double complex :: MomTen(5), CC, Pmu
274  double precision :: MomMax
275  integer :: IndsCoef(2),mu1,mu2,nsum,mu,nu,i,a,cnt,r0
276 
277 
278  tc = 0d0
279  tc(1) = coefsc(0,0,0)
280  tcerr = 0d0
281  tcerr(0) = coefscerr(0)
282 
283  do nsum=1,rmax/2
284  cc = coefsc(nsum,0,0)
285 
286  do nu=rts(nsum-1)+1,rts(nsum)
287  tc(addgtab(1,nu)) = tc(addgtab(1,nu)) + cc*cftab(1,nu)
288  end do
289 
290  end do
291 
292  tcuv = 0d0
293  if (calcuv_ten) then
294 ! TCuv(1) = CoefsC(0,0,0)
295 
296  do nsum=1,rmax/2
297  cc = coefscuv(nsum,0,0)
298 
299  do nu=rts(nsum-1)+1,rts(nsum)
300  tcuv(addgtab(1,nu)) = tcuv(addgtab(1,nu)) + cc*cftab(1,nu)
301  end do
302 
303  end do
304  end if
305 
306 
307  if (rmax.gt.0) then
308 
309  indscoef = (/ 1,0 /)
310  momten(2:5) = momvec(0:3,1)
311  call addtotensorc(1,momten,indscoef)
312 
313  indscoef = (/ 0,1 /)
314  momten(2:5) = momvec(0:3,2)
315  call addtotensorc(1,momten,indscoef)
316 
317  mommax = maxval(abs(momvec))
318  do r0=1,rmax
319  tcerr(r0) = coefscerr(r0)*mommax**r0
320  end do
321 
322  end if
323 
324 
325 
326  contains
327 
328  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
329  ! subroutine AddToTensorC
330  !
331  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
332 
333  recursive subroutine addtotensorc(r,MomTenRec,IndsCoef)
334 
335  integer, intent(in) :: r
336  double complex, intent(in) :: MomTenRec(RtS(r))
337  integer, intent(inout) :: IndsCoef(2)
338  double complex :: MomTen(RtS(r+1)), CC, Pmu
339  integer :: mu1,mu2,nsum,mu,nu,i,a,cnt
340 
341  cc = coefsc(0,indscoef(1),indscoef(2))
342 
343  mu1 = rts(r-1)+1
344  mu2 = rts(r)
345 
346  tc(mu1:mu2) = tc(mu1:mu2) + cc*momtenrec(mu1:mu2)
347 
348  do nsum=1,(rmax-r)/2
349  cc = coefsc(nsum,indscoef(1),indscoef(2))
350 
351  do nu=rts(nsum-1)+1,rts(nsum)
352  do mu=mu1,mu2
353  tc(addgtab(mu,nu)) = tc(addgtab(mu,nu)) + cc*momtenrec(mu)*cftab(mu,nu)
354  end do
355  end do
356 
357  end do
358 
359  if (calcuv_ten) then
360 ! CC = CoefsCuv(0,IndsCoef(1),IndsCoef(2))
361 
362 ! TCuv(mu1:mu2) = TCuv(mu1:mu2) + CC*MomTenRec(mu1:mu2)
363 
364  do nsum=1,(rmax-r)/2
365  cc = coefscuv(nsum,indscoef(1),indscoef(2))
366 
367  do nu=rts(nsum-1)+1,rts(nsum)
368  do mu=mu1,mu2
369  tcuv(addgtab(mu,nu)) = tcuv(addgtab(mu,nu)) + cc*momtenrec(mu)*cftab(mu,nu)
370  end do
371  end do
372 
373  end do
374  end if
375 
376 
377  if (r.lt.rmax) then
378 
379  do i=1,2
380  indscoef(i) = indscoef(i)+1
381 
382  cnt = mu2+1
383  do mu=0,3
384  pmu = momvec(mu,i)
385  do a = mu2-binomtable(r,r+3-mu)+1,mu2
386  momten(cnt)=momtenrec(a)*pmu
387  cnt = cnt+1
388  end do
389  end do
390 
391  call addtotensorc(r+1,momten,indscoef)
392 
393  indscoef(i) = indscoef(i)-1
394  end do
395 
396  end if
397 
398  end subroutine addtotensorc
399 
400 

◆ calctensord()

subroutine buildtensors::calctensord ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TD,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TDuv,
double precision, dimension(0:rmax), intent(out)  TDerr,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax), intent(in)  CoefsD,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax), intent(in)  CoefsDuv,
double precision, dimension(0:rmax), intent(in)  CoefsDerr,
double complex, dimension(0:3,3), intent(in)  MomVec,
integer, intent(in)  rmax 
)

Definition at line 605 of file BuildTensors.F90.

605 
606  integer, intent(in) :: rmax
607  double complex, intent(in) :: MomVec(0:3,3)
608  double complex, intent(in) :: CoefsD(0:rmax/2,0:rmax,0:rmax,0:rmax)
609  double complex, intent(in) :: CoefsDuv(0:rmax/2,0:rmax,0:rmax,0:rmax)
610  double precision, intent(in) :: CoefsDerr(0:rmax)
611  double complex, intent(out) :: TD(0:rmax,0:rmax,0:rmax,0:rmax)
612  double complex, intent(out) :: TDuv(0:rmax,0:rmax,0:rmax,0:rmax)
613  double precision, intent(out) :: TDerr(0:rmax)
614  double complex :: TD_aux(RtS(rmax)), TDuv_aux(RtS(rmax))
615  integer :: mu
616 
617  call calctensord_list(td_aux,tduv_aux,tderr,coefsd,coefsduv,coefsderr,momvec,rmax)
618 
619  do mu=1,rts(rmax)
620  td(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = td_aux(mu)
621  end do
622 
623  if (calcuv_ten) then
624  do mu=1,rts(rmax)
625  tduv(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = tduv_aux(mu)
626  end do
627  else
628  tduv = 0d0
629  end if
630 

◆ calctensord_list()

subroutine buildtensors::calctensord_list ( double complex, dimension(rts(rmax)), intent(out)  TD,
double complex, dimension(rts(rmax)), intent(out)  TDuv,
double precision, dimension(0:rmax), intent(out)  TDerr,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax), intent(in)  CoefsD,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax), intent(in)  CoefsDuv,
double precision, dimension(0:rmax), intent(in)  CoefsDerr,
double complex, dimension(0:3,3), intent(in)  MomVec,
integer, intent(in)  rmax 
)

Definition at line 452 of file BuildTensors.F90.

452 
453  integer, intent(in) :: rmax
454  double complex, intent(in) :: MomVec(0:3,3)
455  double complex, intent(in) :: CoefsD(0:rmax/2,0:rmax,0:rmax,0:rmax)
456  double complex, intent(in) :: CoefsDuv(0:rmax/2,0:rmax,0:rmax,0:rmax)
457  double precision, intent(in) :: CoefsDerr(0:rmax)
458  double complex, intent(out) :: TD(RtS(rmax)), TDuv(RtS(rmax))
459  double precision, intent(out) :: TDerr(0:rmax)
460  double complex :: MomTen(5), CD, Pmu
461  double precision :: MomMax
462  integer :: IndsCoef(3),mu1,mu2,nsum,mu,nu,i,a,cnt,r0
463 
464  td = 0d0
465  td(1) = coefsd(0,0,0,0)
466  tderr = 0d0
467  tderr(0) = coefsderr(0)
468 
469  do nsum=1,rmax/2
470  cd = coefsd(nsum,0,0,0)
471 
472  do nu=rts(nsum-1)+1,rts(nsum)
473  td(addgtab(1,nu)) = td(addgtab(1,nu)) + cd*cftab(1,nu)
474  end do
475 
476  end do
477 
478  tduv = 0d0
479  if (calcuv_ten) then
480 ! TDuv(1) = CoefsDuv(0,0,0,0)
481 
482 ! do nsum=1,rmax/2
483  do nsum=2,rmax/2
484  cd = coefsduv(nsum,0,0,0)
485 
486  do nu=rts(nsum-1)+1,rts(nsum)
487  tduv(addgtab(1,nu)) = tduv(addgtab(1,nu)) + cd*cftab(1,nu)
488  end do
489 
490  end do
491  end if
492 
493 
494  if (rmax.gt.0) then
495 
496  indscoef = (/ 1,0,0 /)
497  momten(2:5) = momvec(0:3,1)
498  call addtotensord(momten,indscoef,1)
499 
500  indscoef = (/ 0,1,0 /)
501  momten(2:5) = momvec(0:3,2)
502  call addtotensord(momten,indscoef,1)
503 
504  indscoef = (/ 0,0,1 /)
505  momten(2:5) = momvec(0:3,3)
506  call addtotensord(momten,indscoef,1)
507 
508  mommax = maxval(abs(momvec))
509  do r0=1,rmax
510  tderr(r0) = coefsderr(r0)*mommax**r0
511  end do
512 
513  end if
514 
515 
516 
517  contains
518 
519  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
520  ! subroutine AddToTensorD
521  !
522  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
523 
524  recursive subroutine addtotensord(MomTenRec,IndsCoef,r)
525 
526  integer, intent(in) :: r
527  integer, intent(inout) :: IndsCoef(3)
528  double complex, intent(in) :: MomTenRec(RtS(r))
529  double complex :: MomTen(RtS(r+1)), CD, Pmu
530  integer :: mu1,mu2,nsum,mu,nu,i,a,cnt
531 
532  cd = coefsd(0,indscoef(1),indscoef(2),indscoef(3))
533 
534  mu1 = rts(r-1)+1
535  mu2 = rts(r)
536 
537  td(mu1:mu2) = td(mu1:mu2) + cd*momtenrec(mu1:mu2)
538 
539  do nsum=1,(rmax-r)/2
540  cd = coefsd(nsum,indscoef(1),indscoef(2),indscoef(3))
541 
542  do nu=rts(nsum-1)+1,rts(nsum)
543  do mu=mu1,mu2
544  td(addgtab(mu,nu)) = td(addgtab(mu,nu)) + cd*momtenrec(mu)*cftab(mu,nu)
545  end do
546  end do
547 
548  end do
549 
550  if (calcuv_ten) then
551 ! CD = CoefsDuv(0,IndsCoef(1),IndsCoef(2),IndsCoef(3))
552 
553 ! TDuv(mu1:mu2) = TDuv(mu1:mu2) + CD*MomTenRec(mu1:mu2)
554 
555 ! do nsum=1,(rmax-r)/2
556  do nsum=2,(rmax-r)/2
557  cd = coefsduv(nsum,indscoef(1),indscoef(2),indscoef(3))
558 
559  do nu=rts(nsum-1)+1,rts(nsum)
560  do mu=mu1,mu2
561  tduv(addgtab(mu,nu)) = tduv(addgtab(mu,nu)) &
562  + cd*momtenrec(mu)*cftab(mu,nu)
563  end do
564  end do
565 
566  end do
567  end if
568 
569 
570  if (r.lt.rmax) then
571 
572  do i=1,3
573  indscoef(i) = indscoef(i)+1
574 
575  cnt = mu2+1
576  do mu=0,3
577  pmu = momvec(mu,i)
578  do a = mu2-binomtable(r,r+3-mu)+1,mu2
579  momten(cnt)=momtenrec(a)*pmu
580  cnt = cnt+1
581  end do
582  end do
583 
584  call addtotensord(momten,indscoef,r+1)
585 
586  indscoef(i) = indscoef(i)-1
587  end do
588 
589  end if
590 
591  end subroutine addtotensord
592 

◆ calctensore()

subroutine buildtensors::calctensore ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TE,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TEuv,
double precision, dimension(0:rmax), intent(out)  TEerr,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsE,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsEuv,
double precision, dimension(0:rmax), intent(in)  CoefsEerr,
double complex, dimension(0:3,4), intent(in)  MomVec,
integer, intent(in)  rmax 
)

Definition at line 804 of file BuildTensors.F90.

804 
805  integer, intent(in) :: rmax
806  double complex, intent(in) :: MomVec(0:3,4)
807  double complex, intent(in) :: CoefsE(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
808  double complex, intent(in) :: CoefsEuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
809  double precision, intent(in) :: CoefsEerr(0:rmax)
810  double complex, intent(out) :: TE(0:rmax,0:rmax,0:rmax,0:rmax)
811  double complex, intent(out) :: TEuv(0:rmax,0:rmax,0:rmax,0:rmax)
812  double precision, intent(out) :: TEerr(0:rmax)
813  double complex :: TE_aux(RtS(rmax)), TEuv_aux(RtS(rmax))
814  integer :: mu
815 
816  call calctensore_list(te_aux,teuv_aux,teerr,coefse,coefseuv,coefseerr,momvec,rmax)
817 
818  do mu=1,rts(rmax)
819  te(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = te_aux(mu)
820  end do
821 
822  if (calcuv_ten) then
823  do mu=1,rts(rmax)
824  teuv(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = teuv_aux(mu)
825  end do
826  else
827  teuv = 0d0
828  end if
829 
830 

◆ calctensore_list()

subroutine buildtensors::calctensore_list ( double complex, dimension(rts(rmax)), intent(out)  TE,
double complex, dimension(rts(rmax)), intent(out)  TEuv,
double precision, dimension(0:rmax), intent(out)  TEerr,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsE,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsEuv,
double precision, dimension(0:rmax), intent(in)  CoefsEerr,
double complex, dimension(0:3,4), intent(in)  MomVec,
integer, intent(in)  rmax 
)

Definition at line 643 of file BuildTensors.F90.

643 
644  integer, intent(in) :: rmax
645  double complex, intent(in) :: MomVec(0:3,4)
646  double complex, intent(in) :: CoefsE(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
647  double complex, intent(in) :: CoefsEuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
648  double precision, intent(in) :: CoefsEerr(0:rmax)
649  double complex, intent(out) :: TE(RtS(rmax)), TEuv(RtS(rmax))
650  double precision, intent(out) :: TEerr(0:rmax)
651  double complex :: MomTen(5), CE, Pmu
652  double precision :: MomMax
653  integer :: IndsCoef(4),mu1,mu2,nsum,mu,nu,i,a,cnt,r0
654 
655 
656  te = 0d0
657  te(1) = coefse(0,0,0,0,0)
658  teerr = 0d0
659  teerr(0) = coefseerr(0)
660 
661  do nsum=1,rmax/2
662  ce = coefse(nsum,0,0,0,0)
663 
664  do nu=rts(nsum-1)+1,rts(nsum)
665  te(addgtab(1,nu)) = te(addgtab(1,nu)) + ce*cftab(1,nu)
666  end do
667 
668  end do
669 
670 
671  teuv = 0d0
672  if (calcuv_ten) then
673 ! TEuv(1) = CoefsEuv(0,0,0,0,0)
674 
675 ! do nsum=1,rmax/2
676  do nsum=3,rmax/2
677  ce = coefseuv(nsum,0,0,0,0)
678 
679  do nu=rts(nsum-1)+1,rts(nsum)
680  teuv(addgtab(1,nu)) = teuv(addgtab(1,nu)) + ce*cftab(1,nu)
681  end do
682 
683  end do
684  end if
685 
686 
687  if (rmax.gt.0) then
688 
689  indscoef = (/ 1,0,0,0 /)
690  momten(2:5) = momvec(0:3,1)
691  call addtotensore(momten,indscoef,1)
692 
693  indscoef = (/ 0,1,0,0 /)
694  momten(2:5) = momvec(0:3,2)
695  call addtotensore(momten,indscoef,1)
696 
697  indscoef = (/ 0,0,1,0 /)
698  momten(2:5) = momvec(0:3,3)
699  call addtotensore(momten,indscoef,1)
700 
701  indscoef = (/ 0,0,0,1 /)
702  momten(2:5) = momvec(0:3,4)
703  call addtotensore(momten,indscoef,1)
704 
705  mommax = maxval(abs(momvec))
706  do r0=1,rmax
707  teerr(r0) = coefseerr(r0)*mommax**r0
708  end do
709 
710  end if
711 
712 
713 
714  contains
715 
716  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
717  ! subroutine AddToTensorE
718  !
719  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
720 
721  recursive subroutine addtotensore(MomTenRec,IndsCoef,r)
722 
723  integer, intent(in) :: r
724  integer, intent(inout) :: IndsCoef(4)
725  double complex, intent(in) :: MomTenRec(RtS(r))
726  double complex :: MomTen(RtS(r+1)), CE, Pmu
727  double precision :: CEmax
728  integer :: mu1,mu2,nsum,mu,nu,i,a,cnt
729 
730  ce = coefse(0,indscoef(1),indscoef(2),indscoef(3),indscoef(4))
731  cemax = abs(ce) ! abs(CoefsE(0,0,0,0,0))
732 
733  mu1 = rts(r-1)+1
734  mu2 = rts(r)
735 
736  te(mu1:mu2) = te(mu1:mu2) + ce*momtenrec(mu1:mu2)
737 
738  do nsum=1,(rmax-r)/2
739  ce = coefse(nsum,indscoef(1),indscoef(2),indscoef(3),indscoef(4))
740 
741  do nu=rts(nsum-1)+1,rts(nsum)
742  do mu=mu1,mu2
743  te(addgtab(mu,nu)) = te(addgtab(mu,nu)) + ce*momtenrec(mu)*cftab(mu,nu)
744  end do
745  end do
746 
747  end do
748 
749  if (calcuv_ten) then
750 ! CE = CoefsEuv(0,IndsCoef(1),IndsCoef(2),IndsCoef(3),IndsCoef(4))
751 
752 ! TEuv(mu1:mu2) = TEuv(mu1:mu2) + CE*MomTenRec(mu1:mu2)
753 
754 ! do nsum=1,(rmax-r)/2
755  do nsum=1,(rmax-r)/2
756  ce = coefseuv(nsum,indscoef(1),indscoef(2),indscoef(3),indscoef(4))
757 
758  do nu=rts(nsum-1)+1,rts(nsum)
759  do mu=mu1,mu2
760  teuv(addgtab(mu,nu)) = teuv(addgtab(mu,nu)) &
761  + ce*momtenrec(mu)*cftab(mu,nu)
762  end do
763  end do
764 
765  end do
766  end if
767 
768 
769  if (r.lt.rmax) then
770 
771  do i=1,4
772  indscoef(i) = indscoef(i)+1
773 
774  cnt = mu2+1
775  do mu=0,3
776  pmu = momvec(mu,i)
777  do a = mu2-binomtable(r,r+3-mu)+1,mu2
778  momten(cnt)=momtenrec(a)*pmu
779  cnt = cnt+1
780  end do
781  end do
782 
783  call addtotensore(momten,indscoef,r+1)
784 
785  indscoef(i) = indscoef(i)-1
786  end do
787 
788  end if
789 
790  end subroutine addtotensore
791 

◆ calctensorf()

subroutine buildtensors::calctensorf ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TF,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TFuv,
double precision, dimension(0:rmax), intent(out)  TFerr,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsF,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsFuv,
double precision, dimension(0:rmax), intent(in)  CoefsFerr,
double complex, dimension(0:3,5), intent(in)  MomVec,
integer, intent(in)  rmax 
)

Definition at line 1004 of file BuildTensors.F90.

1004 
1005  integer, intent(in) :: rmax
1006  double complex, intent(in) :: MomVec(0:3,5)
1007  double complex, intent(in) :: CoefsF(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
1008  double complex, intent(in) :: CoefsFuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
1009  double precision, intent(in) :: CoefsFerr(0:rmax)
1010  double complex, intent(out) :: TF(0:rmax,0:rmax,0:rmax,0:rmax)
1011  double complex, intent(out) :: TFuv(0:rmax,0:rmax,0:rmax,0:rmax)
1012  double precision, intent(out) :: TFerr(0:rmax)
1013  double complex :: TF_aux(RtS(rmax)), TFuv_aux(RtS(rmax))
1014  integer :: mu
1015 
1016  call calctensorf_list(tf_aux,tfuv_aux,tferr,coefsf,coefsfuv,coefsferr,momvec,rmax)
1017 
1018  do mu=1,rts(rmax)
1019  tf(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = tf_aux(mu)
1020  end do
1021 
1022  if (calcuv_ten) then
1023  do mu=1,rts(rmax)
1024  tfuv(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = tfuv_aux(mu)
1025  end do
1026  else
1027  tfuv = 0d0
1028  end if
1029 

◆ calctensorf_list()

subroutine buildtensors::calctensorf_list ( double complex, dimension(rts(rmax)), intent(out)  TF,
double complex, dimension(rts(rmax)), intent(out)  TFuv,
double precision, dimension(0:rmax), intent(out)  TFerr,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsF,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsFuv,
double precision, dimension(0:rmax), intent(in)  CoefsFerr,
double complex, dimension(0:3,5), intent(in)  MomVec,
integer, intent(in)  rmax 
)

Definition at line 843 of file BuildTensors.F90.

843 
844  integer, intent(in) :: rmax
845  double complex, intent(in) :: MomVec(0:3,5)
846  double complex, intent(in) :: CoefsF(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
847  double complex, intent(in) :: CoefsFuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
848  double precision, intent(in) :: CoefsFerr(0:rmax)
849  double complex, intent(out) :: TF(RtS(rmax)), TFuv(RtS(rmax))
850  double precision, intent(out) :: TFerr(0:rmax)
851  double precision :: MomMax
852  double complex :: MomTen(5), CF, Pmu
853  integer :: IndsCoef(5),mu1,mu2,nsum,mu,nu,i,a,cnt,r0
854 
855  tf = 0d0
856  tf(1) = coefsf(0,0,0,0,0,0)
857  tferr = 0d0
858  tferr(0) = coefsferr(0)
859 
860  do nsum=1,rmax/2
861  cf = coefsf(nsum,0,0,0,0,0)
862 
863  do nu=rts(nsum-1)+1,rts(nsum)
864  tf(addgtab(1,nu)) = tf(addgtab(1,nu)) + cf*cftab(1,nu)
865  end do
866 
867  end do
868 
869  tfuv = 0d0
870  if (calcuv_ten) then
871 ! TFuv(1) = CoefsFuv(0,0,0,0,0,0)
872 
873 ! do nsum=1,rmax/2
874  do nsum=4,rmax/2
875  cf = coefsfuv(nsum,0,0,0,0,0)
876 
877  do nu=rts(nsum-1)+1,rts(nsum)
878  tfuv(addgtab(1,nu)) = tfuv(addgtab(1,nu)) + cf*cftab(1,nu)
879  end do
880 
881  end do
882  end if
883 
884 
885  if (rmax.gt.0) then
886 
887  indscoef = (/ 1,0,0,0,0 /)
888  momten(2:5) = momvec(0:3,1)
889  call addtotensorf(momten,indscoef,1)
890 
891  indscoef = (/ 0,1,0,0,0 /)
892  momten(2:5) = momvec(0:3,2)
893  call addtotensorf(momten,indscoef,1)
894 
895  indscoef = (/ 0,0,1,0,0 /)
896  momten(2:5) = momvec(0:3,3)
897  call addtotensorf(momten,indscoef,1)
898 
899  indscoef = (/ 0,0,0,1,0 /)
900  momten(2:5) = momvec(0:3,4)
901  call addtotensorf(momten,indscoef,1)
902 
903  indscoef = (/ 0,0,0,0,1 /)
904  momten(2:5) = momvec(0:3,5)
905  call addtotensorf(momten,indscoef,1)
906 
907  mommax = maxval(abs(momvec))
908  do r0=1,rmax
909  tferr(r0) = coefsferr(r0)*mommax**r0
910  end do
911 
912  end if
913 
914 
915 
916  contains
917 
918  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
919  ! subroutine AddToTensorF
920  !
921  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
922 
923  recursive subroutine addtotensorf(MomTenRec,IndsCoef,r)
924 
925  integer, intent(in) :: r
926  integer, intent(inout) :: IndsCoef(5)
927  double complex, intent(in) :: MomTenRec(RtS(r))
928  double complex :: MomTen(RtS(r+1)), CF, Pmu
929  integer :: mu1,mu2,nsum,mu,nu,i,a,cnt
930 
931  cf = coefsf(0,indscoef(1),indscoef(2),indscoef(3),indscoef(4),indscoef(5))
932 
933  mu1 = rts(r-1)+1
934  mu2 = rts(r)
935 
936  tf(mu1:mu2) = tf(mu1:mu2) + cf*momtenrec(mu1:mu2)
937 
938  do nsum=1,(rmax-r)/2
939  cf = coefsf(nsum,indscoef(1),indscoef(2),indscoef(3),indscoef(4),indscoef(5))
940 
941  do nu=rts(nsum-1)+1,rts(nsum)
942  do mu=mu1,mu2
943  tf(addgtab(mu,nu)) = tf(addgtab(mu,nu)) + cf*momtenrec(mu)*cftab(mu,nu)
944  end do
945  end do
946 
947  end do
948 
949  if (calcuv_ten) then
950 ! CF = CoefsFuv(0,IndsCoef(1),IndsCoef(2),IndsCoef(3),IndsCoef(4),IndsCoef(5))
951 
952 ! TFuv(mu1:mu2) = TFuv(mu1:mu2) + CF*MomTenRec(mu1:mu2)
953 
954 ! do nsum=1,(rmax-r)/2
955  do nsum=4,(rmax-r)/2
956  cf = coefsfuv(nsum,indscoef(1),indscoef(2),indscoef(3),indscoef(4),indscoef(5))
957 
958  do nu=rts(nsum-1)+1,rts(nsum)
959  do mu=mu1,mu2
960  tfuv(addgtab(mu,nu)) = tfuv(addgtab(mu,nu)) &
961  + cf*momtenrec(mu)*cftab(mu,nu)
962  end do
963  end do
964 
965  end do
966  end if
967 
968 
969  if (r.lt.rmax) then
970 
971  do i=1,5
972  indscoef(i) = indscoef(i)+1
973 
974  cnt = mu2+1
975  do mu=0,3
976  pmu = momvec(mu,i)
977  do a = mu2-binomtable(r,r+3-mu)+1,mu2
978  momten(cnt)=momtenrec(a)*pmu
979  cnt = cnt+1
980  end do
981  end do
982 
983  call addtotensorf(momten,indscoef,r+1)
984 
985  indscoef(i) = indscoef(i)-1
986  end do
987 
988  end if
989 
990  end subroutine addtotensorf
991 

◆ calctensorfuv()

subroutine buildtensors::calctensorfuv ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TFuv,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsFuv,
double complex, dimension(0:3,5), intent(in)  MomVec,
integer, intent(in)  rmax 
)

Definition at line 1153 of file BuildTensors.F90.

1153 
1154  integer, intent(in) :: rmax
1155  double complex, intent(in) :: MomVec(0:3,5)
1156  double complex, intent(in) :: CoefsFuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
1157  double complex, intent(out) :: TFuv(0:rmax,0:rmax,0:rmax,0:rmax)
1158  double complex :: TFuv_aux(RtS(rmax))
1159  integer :: mu
1160 
1161  call calctensorfuv_list(tfuv_aux,coefsfuv,momvec,rmax)
1162 
1163  do mu=1,rts(rmax)
1164  tfuv(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = tfuv_aux(mu)
1165  end do
1166 

◆ calctensorfuv_list()

subroutine buildtensors::calctensorfuv_list ( double complex, dimension(rts(rmax)), intent(out)  TFuv,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsFuv,
double complex, dimension(0:3,5), intent(in)  MomVec,
integer, intent(in)  rmax 
)

Definition at line 1042 of file BuildTensors.F90.

1042 
1043  integer, intent(in) :: rmax
1044  double complex, intent(in) :: MomVec(0:3,5)
1045  double complex, intent(in) :: CoefsFuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
1046  double complex, intent(out) :: TFuv(RtS(rmax))
1047  double complex :: MomTen(5), CF, Pmu
1048  integer :: IndsCoef(5),mu1,mu2,nsum,mu,nu,i,a,cnt
1049 
1050 
1051  tfuv = 0d0
1052  do nsum=4,rmax/2
1053  cf = coefsfuv(nsum,0,0,0,0,0)
1054 
1055  do nu=rts(nsum-1)+1,rts(nsum)
1056  tfuv(addgtab(1,nu)) = tfuv(addgtab(1,nu)) + cf*cftab(1,nu)
1057  end do
1058 
1059  end do
1060 
1061  if (rmax.gt.8) then
1062 
1063  indscoef = (/ 1,0,0,0,0 /)
1064  momten(2:5) = momvec(0:3,1)
1065  call addtotensorfuv(momten,indscoef,1)
1066 
1067  indscoef = (/ 0,1,0,0,0 /)
1068  momten(2:5) = momvec(0:3,2)
1069  call addtotensorfuv(momten,indscoef,1)
1070 
1071  indscoef = (/ 0,0,1,0,0 /)
1072  momten(2:5) = momvec(0:3,3)
1073  call addtotensorfuv(momten,indscoef,1)
1074 
1075  indscoef = (/ 0,0,0,1,0 /)
1076  momten(2:5) = momvec(0:3,4)
1077  call addtotensorfuv(momten,indscoef,1)
1078 
1079  indscoef = (/ 0,0,0,0,1 /)
1080  momten(2:5) = momvec(0:3,5)
1081  call addtotensorfuv(momten,indscoef,1)
1082 
1083  end if
1084 
1085 
1086 
1087  contains
1088 
1089  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1090  ! subroutine AddToTensorFuv
1091  !
1092  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1093 
1094  recursive subroutine addtotensorfuv(MomTenRec,IndsCoef,r)
1095 
1096  integer, intent(in) :: r
1097  integer, intent(inout) :: IndsCoef(5)
1098  double complex, intent(in) :: MomTenRec(RtS(r))
1099  double complex :: MomTen(RtS(r+1)), CF, Pmu
1100  integer :: mu1,mu2,nsum,mu,nu,i,a,cnt
1101 
1102  mu1 = rts(r-1)+1
1103  mu2 = rts(r)
1104 
1105  do nsum=4,(rmax-r)/2
1106  cf = coefsfuv(nsum,indscoef(1),indscoef(2),indscoef(3),indscoef(4),indscoef(5))
1107 
1108  do nu=rts(nsum-1)+1,rts(nsum)
1109  do mu=mu1,mu2
1110  tfuv(addgtab(mu,nu)) = tfuv(addgtab(mu,nu)) &
1111  + cf*momtenrec(mu)*cftab(mu,nu)
1112  end do
1113  end do
1114 
1115  end do
1116 
1117 
1118  if (r.lt.rmax-8) then
1119 
1120  do i=1,5
1121  indscoef(i) = indscoef(i)+1
1122 
1123  cnt = mu2+1
1124  do mu=0,3
1125  pmu = momvec(mu,i)
1126  do a = mu2-binomtable(r,r+3-mu)+1,mu2
1127  momten(cnt)=momtenrec(a)*pmu
1128  cnt = cnt+1
1129  end do
1130  end do
1131 
1132  call addtotensorfuv(momten,indscoef,r+1)
1133 
1134  indscoef(i) = indscoef(i)-1
1135  end do
1136 
1137  end if
1138 
1139  end subroutine addtotensorfuv
1140 

◆ calctensorg()

subroutine buildtensors::calctensorg ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TG,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TGuv,
double precision, dimension(0:rmax), intent(out)  TGerr,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsG,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsGuv,
double precision, dimension(0:rmax), intent(in)  CoefsGerr,
double complex, dimension(0:3,6), intent(in)  MomVec,
integer, intent(in)  rmax 
)

Definition at line 1344 of file BuildTensors.F90.

1344 
1345  integer, intent(in) :: rmax
1346  double complex, intent(in) :: MomVec(0:3,6)
1347  double complex, intent(in) :: CoefsG(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
1348  double complex, intent(in) :: CoefsGuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
1349  double precision, intent(in) :: CoefsGerr(0:rmax)
1350  double complex, intent(out) :: TG(0:rmax,0:rmax,0:rmax,0:rmax)
1351  double complex, intent(out) :: TGuv(0:rmax,0:rmax,0:rmax,0:rmax)
1352  double precision, intent(out) :: TGerr(0:rmax)
1353  double complex :: TG_aux(RtS(rmax)), TGuv_aux(RtS(rmax))
1354  integer :: mu
1355 
1356  call calctensorg_list(tg_aux,tguv_aux,tgerr,coefsg,coefsguv,coefsgerr,momvec,rmax)
1357 
1358  do mu=1,rts(rmax)
1359  tg(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = tg_aux(mu)
1360  end do
1361 
1362  if (calcuv_ten) then
1363  do mu=1,rts(rmax)
1364  tguv(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = tguv_aux(mu)
1365  end do
1366  else
1367  tguv = 0d0
1368  end if
1369 
1370 

◆ calctensorg_list()

subroutine buildtensors::calctensorg_list ( double complex, dimension(rts(rmax)), intent(out)  TG,
double complex, dimension(rts(rmax)), intent(out)  TGuv,
double precision, dimension(0:rmax), intent(out)  TGerr,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsG,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(in)  CoefsGuv,
double precision, dimension(0:rmax), intent(in)  CoefsGerr,
double complex, dimension(0:3,6), intent(in)  MomVec,
integer, intent(in)  rmax 
)

Definition at line 1178 of file BuildTensors.F90.

1178 
1179  integer, intent(in) :: rmax
1180  double complex, intent(in) :: MomVec(0:3,6)
1181  double complex, intent(in) :: CoefsG(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
1182  double complex, intent(in) :: CoefsGuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
1183  double precision, intent(in) :: CoefsGerr(0:rmax)
1184  double complex, intent(out) :: TG(RtS(rmax)), TGuv(RtS(rmax))
1185  double precision, intent(out) :: TGerr(0:rmax)
1186  double complex :: MomTen(5), CG, Pmu
1187  double precision :: MomMax
1188  integer :: IndsCoef(6),mu1,mu2,nsum,mu,nu,i,a,cnt,r0
1189 
1190 
1191  tg = 0d0
1192  tg(1) = coefsg(0,0,0,0,0,0,0)
1193  tgerr = 0d0
1194  tgerr(0) = coefsgerr(0)
1195 
1196  do nsum=1,rmax/2
1197  cg = coefsg(nsum,0,0,0,0,0,0)
1198 
1199  do nu=rts(nsum-1)+1,rts(nsum)
1200  tg(addgtab(1,nu)) = tg(addgtab(1,nu)) + cg*cftab(1,nu)
1201  end do
1202 
1203  end do
1204 
1205  tguv = 0d0
1206  if (calcuv_ten) then
1207 ! TGuv(1) = CoefsGuv(0,0,0,0,0,0,0)
1208 
1209 ! do nsum=1,rmax/2
1210  do nsum=5,rmax/2
1211  cg = coefsguv(nsum,0,0,0,0,0,0)
1212 
1213  do nu=rts(nsum-1)+1,rts(nsum)
1214  tguv(addgtab(1,nu)) = tguv(addgtab(1,nu)) + cg*cftab(1,nu)
1215  end do
1216 
1217  end do
1218  end if
1219 
1220 
1221  if (rmax.gt.0) then
1222 
1223  indscoef = (/ 1,0,0,0,0,0 /)
1224  momten(2:5) = momvec(0:3,1)
1225  call addtotensorg(momten,indscoef,1)
1226 
1227  indscoef = (/ 0,1,0,0,0,0 /)
1228  momten(2:5) = momvec(0:3,2)
1229  call addtotensorg(momten,indscoef,1)
1230 
1231  indscoef = (/ 0,0,1,0,0,0 /)
1232  momten(2:5) = momvec(0:3,3)
1233  call addtotensorg(momten,indscoef,1)
1234 
1235  indscoef = (/ 0,0,0,1,0,0 /)
1236  momten(2:5) = momvec(0:3,4)
1237  call addtotensorg(momten,indscoef,1)
1238 
1239  indscoef = (/ 0,0,0,0,1,0 /)
1240  momten(2:5) = momvec(0:3,5)
1241  call addtotensorg(momten,indscoef,1)
1242 
1243  indscoef = (/ 0,0,0,0,0,1 /)
1244  momten(2:5) = momvec(0:3,6)
1245  call addtotensorg(momten,indscoef,1)
1246 
1247  mommax = maxval(abs(momvec))
1248  do r0=1,rmax
1249  tgerr(r0) = coefsgerr(r0)*mommax**r0
1250  end do
1251 
1252  end if
1253 
1254 
1255 
1256  contains
1257 
1258  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1259  ! subroutine AddToTensorG
1260  !
1261  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1262 
1263  recursive subroutine addtotensorg(MomTenRec,IndsCoef,r)
1264 
1265  integer :: r
1266  integer, intent(inout) :: IndsCoef(6)
1267  double complex, intent(in) :: MomTenRec(RtS(r))
1268  double complex :: MomTen(RtS(r+1)), CG, Pmu
1269  integer :: mu1,mu2,nsum,mu,nu,i,a,cnt
1270 
1271  cg = coefsg(0,indscoef(1),indscoef(2),indscoef(3),indscoef(4),indscoef(5),indscoef(6))
1272 
1273  mu1 = rts(r-1)+1
1274  mu2 = rts(r)
1275 
1276  tg(mu1:mu2) = tg(mu1:mu2) + cg*momtenrec(mu1:mu2)
1277 
1278  do nsum=1,(rmax-r)/2
1279  cg = coefsg(nsum,indscoef(1),indscoef(2),indscoef(3),indscoef(4),indscoef(5),indscoef(6))
1280 
1281  do nu=rts(nsum-1)+1,rts(nsum)
1282  do mu=mu1,mu2
1283  tg(addgtab(mu,nu)) = tg(addgtab(mu,nu)) + cg*momtenrec(mu)*cftab(mu,nu)
1284  end do
1285  end do
1286 
1287  end do
1288 
1289  if (calcuv_ten) then
1290 ! CG = CoefsGuv(0,IndsCoef(1),IndsCoef(2),IndsCoef(3),IndsCoef(4),IndsCoef(5),IndsCoef(6))
1291 
1292 ! TGuv(mu1:mu2) = TGuv(mu1:mu2) + CG*MomTenRec(mu1:mu2)
1293 
1294 ! do nsum=1,(rmax-r)/2
1295  do nsum=5,(rmax-r)/2
1296  cg = coefsguv(nsum,indscoef(1),indscoef(2),indscoef(3),indscoef(4),indscoef(5),indscoef(6))
1297 
1298  do nu=rts(nsum-1)+1,rts(nsum)
1299  do mu=mu1,mu2
1300  tguv(addgtab(mu,nu)) = tguv(addgtab(mu,nu)) &
1301  + cg*momtenrec(mu)*cftab(mu,nu)
1302  end do
1303  end do
1304 
1305  end do
1306  end if
1307 
1308 
1309  if (r.lt.rmax) then
1310 
1311  do i=1,6
1312  indscoef(i) = indscoef(i)+1
1313 
1314  cnt = mu2+1
1315  do mu=0,3
1316  pmu = momvec(mu,i)
1317  do a = mu2-binomtable(r,r+3-mu)+1,mu2
1318  momten(cnt)=momtenrec(a)*pmu
1319  cnt = cnt+1
1320  end do
1321  end do
1322 
1323  call addtotensorg(momten,indscoef,r+1)
1324 
1325  indscoef(i) = indscoef(i)-1
1326  end do
1327 
1328  end if
1329 
1330  end subroutine addtotensorg
1331 

◆ calctensortn()

subroutine buildtensors::calctensortn ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TN,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TNuv,
double precision, dimension(0:rmax), intent(out)  TNerr,
double complex, dimension(ncoefs(rmax,n)), intent(in)  CoefsN,
double complex, dimension(ncoefs(rmax,n)), intent(in)  CoefsNuv,
double precision, dimension(0:rmax), intent(in)  CoefsNerr,
double complex, dimension(0:3,n-1), intent(in)  MomVec,
integer, intent(in)  N,
integer, intent(in)  rmax 
)

Definition at line 1546 of file BuildTensors.F90.

1546 
1547  integer, intent(in) :: N,rmax
1548  double complex, intent(in) :: MomVec(0:3,N-1)
1549  double complex, intent(in) :: CoefsN(NCoefs(rmax,N))
1550  double complex, Intent(in) :: CoefsNuv(NCoefs(rmax,N))
1551  double precision, intent(in) :: CoefsNerr(0:rmax)
1552  double complex, intent(out) :: TN(0:rmax,0:rmax,0:rmax,0:rmax)
1553  double complex, intent(out) :: TNuv(0:rmax,0:rmax,0:rmax,0:rmax)
1554  double precision, intent(out) :: TNerr(0:rmax)
1555  double complex :: TN_aux(RtS(rmax)), TNuv_aux(RtS(rmax))
1556  double complex :: MomTen(5), CN, Pmu
1557  integer :: IndsCoef(5),mu1,mu2,nsum,mu,nu,ind,a,cnt
1558 
1559  call calctensortn_list(tn_aux,tnuv_aux,tnerr,coefsn,coefsnuv,coefsnerr,momvec,n,rmax)
1560 
1561  do mu=1,rts(rmax)
1562  tn(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = tn_aux(mu)
1563  end do
1564 
1565  if (calcuv_ten) then
1566  do mu=1,rts(rmax)
1567  tnuv(lorindtab(0,mu),lorindtab(1,mu),lorindtab(2,mu),lorindtab(3,mu)) = tnuv_aux(mu)
1568  end do
1569  else
1570  tnuv = 0d0
1571  end if
1572 
1573 

◆ calctensortn_list()

subroutine buildtensors::calctensortn_list ( double complex, dimension(rts(rmax)), intent(out)  TN,
double complex, dimension(rts(rmax)), intent(out)  TNuv,
double precision, dimension(0:rmax), intent(out)  TNerr,
double complex, dimension(ncoefs(rmax,n)), intent(in)  CoefsN,
double complex, dimension(ncoefs(rmax,n)), intent(in)  CoefsNuv,
double precision, dimension(0:rmax), intent(in)  CoefsNerr,
double complex, dimension(0:3,n-1), intent(in)  MomVec,
integer, intent(in)  N,
integer, intent(in)  rmax 
)

Definition at line 1383 of file BuildTensors.F90.

1383 
1384  integer, intent(in) :: N,rmax
1385  double complex, intent(in) :: MomVec(0:3,N-1)
1386  double complex, intent(in) :: CoefsN(NCoefs(rmax,N)),CoefsNuv(NCoefs(rmax,N))
1387  double precision, intent(in) :: CoefsNerr(0:rmax)
1388  double complex :: CoefsN_aux(BinomTable(rmax,max(N+rmax-2,0)),0:rmax/2,0:rmax)
1389  double complex :: CoefsNuv_aux(BinomTable(rmax,max(N+rmax-2,0)),0:rmax/2,0:rmax)
1390  double complex, intent(out) :: TN(RtS(rmax)), TNuv(RtS(rmax))
1391  double precision, intent(out) :: TNerr(0:rmax)
1392  double precision :: MomMax
1393  double complex :: MomTen(5), CN, Pmu
1394  integer :: IndsCoef(5),mu1,mu2,nsum,mu,nu,ind,a,cnt,r0,n0,i,r
1395 
1396  cnt = 0
1397  do r=0,rmax
1398  do n0=r/2,0,-1
1399  do i=1,binomtable(r-2*n0,max(n+r-2*n0-2,0))
1400 
1401  cnt = cnt+1
1402  coefsn_aux(i,n0,r) = coefsn(cnt)
1403  coefsnuv_aux(i,n0,r) = coefsnuv(cnt)
1404 
1405  end do
1406  end do
1407  end do
1408 
1409 
1410  tn = 0d0
1411  tn(1) = coefsn_aux(1,0,0)
1412  tnerr = 0d0
1413  tnerr(0) = coefsnerr(0)
1414 
1415  do nsum=1,rmax/2
1416  cn = coefsn_aux(1,nsum,2*nsum)
1417 
1418  do nu=rts(nsum-1)+1,rts(nsum)
1419  tn(addgtab(1,nu)) = tn(addgtab(1,nu)) + cn*cftab(1,nu)
1420  end do
1421 
1422  end do
1423 
1424  tnuv = 0d0
1425  if (calcuv_ten) then
1426 ! TNuv(1) = CoefsNuv_aux(1,0,0)
1427 
1428 ! do nsum=1,rmax/2
1429  do nsum=max(n-2,1),rmax/2
1430  cn = coefsnuv_aux(1,nsum,2*nsum)
1431 
1432  do nu=rts(nsum-1)+1,rts(nsum)
1433  tnuv(addgtab(1,nu)) = tnuv(addgtab(1,nu)) + cn*cftab(1,nu)
1434  end do
1435 
1436  end do
1437  end if
1438 
1439 
1440  if (rmax.gt.0) then
1441 
1442  do ind=1,n-1
1443  momten(2:5) = momvec(0:3,ind)
1444  call addtotensortn(momten,ind,1)
1445  end do
1446 
1447  mommax = maxval(abs(momvec))
1448  do r0=1,rmax
1449  tnerr(r0) = coefsnerr(r0)*mommax**r0
1450  end do
1451  end if
1452 
1453 
1454 
1455  contains
1456 
1457  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1458  ! subroutine AddToTensorTN
1459  !
1460  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1461 
1462  recursive subroutine addtotensortn(MomTenRec,ind,r)
1463 
1464  integer, intent(in) :: r,ind
1465  double complex, intent(in) :: MomTenRec(RtS(r))
1466  double complex :: CN, Pmu
1467  double complex, allocatable :: MomTenIt(:)
1468  integer :: mu1,mu2,nsum,mu,nu,i,a,cnt,nind
1469  double precision :: CNmax
1470 
1471  cn = coefsn_aux(ind,0,r)
1472 
1473  mu1 = rts(r-1)+1
1474  mu2 = rts(r)
1475 
1476  tn(mu1:mu2) = tn(mu1:mu2) + cn*momtenrec(mu1:mu2)
1477 
1478  do nsum=1,(rmax-r)/2
1479  cn = coefsn_aux(ind,nsum,r+2*nsum)
1480 
1481  do nu=rts(nsum-1)+1,rts(nsum)
1482  do mu=mu1,mu2
1483  tn(addgtab(mu,nu)) = tn(addgtab(mu,nu)) + cn*momtenrec(mu)*cftab(mu,nu)
1484  end do
1485  end do
1486 
1487  end do
1488 
1489  if (calcuv_ten) then
1490 
1491  if (n.le.2) then
1492  cn = coefsnuv_aux(ind,0,r)
1493  tnuv(mu1:mu2) = tnuv(mu1:mu2) + cn*momtenrec(mu1:mu2)
1494  end if
1495 
1496 ! do nsum=1,(rmax-r)/2
1497  do nsum=max(n-2,1),(rmax-r)/2
1498  cn = coefsnuv_aux(ind,nsum,r+2*nsum)
1499 
1500  do nu=rts(nsum-1)+1,rts(nsum)
1501  do mu=mu1,mu2
1502  tnuv(addgtab(mu,nu)) = tnuv(addgtab(mu,nu)) &
1503  + cn*momtenrec(mu)*cftab(mu,nu)
1504  end do
1505  end do
1506 
1507  end do
1508  end if
1509 
1510 
1511  if (r.lt.rmax) then
1512 
1513  allocate(momtenit(rts(r+1)))
1514  do i=1,n-1
1515  nind = indspiprod(0,i,ind,n-1)
1516 
1517  cnt = mu2+1
1518  do mu=0,3
1519  pmu = momvec(mu,i)
1520  do a = mu2-binomtable(r,r+3-mu)+1,mu2
1521  momtenit(cnt)=momtenrec(a)*pmu
1522  cnt = cnt+1
1523  end do
1524  end do
1525 
1526  call addtotensortn(momtenit,nind,r+1)
1527 
1528  end do
1529 
1530  end if
1531 
1532  end subroutine addtotensortn
1533 

◆ calctensortnuv_list()

subroutine buildtensors::calctensortnuv_list ( double complex, dimension(rts(rmax)), intent(out)  TNuv,
double complex, dimension(binomtable(rmax-2*n+4,max(rmax-n+2,0)),n-2:rmax/2,2*n-4:rmax), intent(in)  CoefsNuv,
double complex, dimension(0:3,n-1), intent(in)  MomVec,
integer, intent(in)  N,
integer, intent(in)  rmax 
)

Definition at line 1586 of file BuildTensors.F90.

1586 
1587  integer, intent(in) :: N,rmax
1588  double complex, intent(in) :: MomVec(0:3,N-1)
1589  double complex, intent(in) :: CoefsNuv(BinomTable(rmax-2*N+4,max(rmax-N+2,0)),N-2:rmax/2,2*N-4:rmax)
1590  double complex, intent(out) :: TNuv(RtS(rmax))
1591  double complex :: MomTen(5), CN, Pmu
1592  integer :: IndsCoef(5),mu1,mu2,nsum,mu,nu,ind,a,cnt,r,n0,i
1593 
1594 
1595  tnuv = 0d0
1596  tnuv(rts(2*n-5)+1) = coefsnuv(1,n-2,2*n-4)
1597 
1598  do nsum=max(n-2,1),rmax/2
1599  cn = coefsnuv(1,nsum,2*nsum)
1600 
1601  do nu=rts(nsum-1)+1,rts(nsum)
1602  tnuv(addgtab(1,nu)) = tnuv(addgtab(1,nu)) + cn*cftab(1,nu)
1603  end do
1604  end do
1605 
1606 
1607  if (rmax.gt.2*n-4) then
1608 
1609  do ind=1,n-1
1610  momten(2:5) = momvec(0:3,ind)
1611  call addtotensortnuv(momten,ind,1)
1612  end do
1613 
1614  end if
1615 
1616 
1617 
1618  contains
1619 
1620  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1621  ! subroutine AddToTensorTNuv
1622  !
1623  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1624 
1625  recursive subroutine addtotensortnuv(MomTenRec,ind,r)
1626 
1627  integer, intent(in) :: r,ind
1628  double complex, intent(in) :: MomTenRec(RtS(r))
1629  double complex :: MomTen(RtS(r+1)), CN, Pmu
1630  integer :: mu1,mu2,nsum,mu,nu,i,a,cnt,nind
1631  double precision :: CNmax
1632 
1633  mu1 = rts(r-1)+1
1634  mu2 = rts(r)
1635 
1636  if (n.le.2) then
1637  cn = coefsnuv(ind,0,r)
1638  tnuv(mu1:mu2) = tnuv(mu1:mu2) + cn*momtenrec(mu1:mu2)
1639  end if
1640 
1641  do nsum=max(n-2,1),(rmax-r)/2
1642  cn = coefsnuv(ind,nsum,r+2*nsum)
1643 
1644  do nu=rts(nsum-1)+1,rts(nsum)
1645  do mu=mu1,mu2
1646  tnuv(addgtab(mu,nu)) = tnuv(addgtab(mu,nu)) &
1647  + cn*momtenrec(mu)*cftab(mu,nu)
1648  end do
1649  end do
1650 
1651  end do
1652 
1653 
1654  if (r.lt.rmax-2*n+4) then
1655 
1656  do i=1,n-1
1657  nind = indspiprod(0,i,ind,n-1)
1658 
1659  cnt = mu2+1
1660  do mu=0,3
1661  pmu = momvec(mu,i)
1662  do a = mu2-binomtable(r,r+3-mu)+1,mu2
1663  momten(cnt)=momtenrec(a)*pmu
1664  cnt = cnt+1
1665  end do
1666  end do
1667 
1668  call addtotensortnuv(momten,nind,r+1)
1669 
1670  end do
1671 
1672  end if
1673 
1674  end subroutine addtotensortnuv
1675 
addtotensorf
recursive subroutine addtotensorf(MomTenRec, IndsCoef, r)
Definition: BuildTensors.F90:924
addtotensore
recursive subroutine addtotensore(MomTenRec, IndsCoef, r)
Definition: BuildTensors.F90:722
addtotensorc
recursive subroutine addtotensorc(r, MomTenRec, IndsCoef)
Definition: BuildTensors.F90:334
addtotensorg
recursive subroutine addtotensorg(MomTenRec, IndsCoef, r)
Definition: BuildTensors.F90:1264
addtotensortnuv
recursive subroutine addtotensortnuv(MomTenRec, ind, r)
Definition: BuildTensors.F90:1626
addtotensord
recursive subroutine addtotensord(MomTenRec, IndsCoef, r)
Definition: BuildTensors.F90:525
addtotensorfuv
recursive subroutine addtotensorfuv(MomTenRec, IndsCoef, r)
Definition: BuildTensors.F90:1095
addtotensortn
recursive subroutine addtotensortn(MomTenRec, ind, r)
Definition: BuildTensors.F90:1463