JHUGen MELA  JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
BuildTensors.F90 File Reference

Go to the source code of this file.

Modules

module  buildtensors
 

Functions/Subroutines

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

Function/Subroutine Documentation

◆ addtotensorc()

recursive subroutine calctensorc_list::addtotensorc ( integer, intent(in)  r,
double complex, dimension(rts(r)), intent(in)  MomTenRec,
integer, dimension(2), intent(inout)  IndsCoef 
)

Definition at line 334 of file BuildTensors.F90.

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 

◆ addtotensord()

recursive subroutine calctensord_list::addtotensord ( double complex, dimension(rts(r)), intent(in)  MomTenRec,
integer, dimension(3), intent(inout)  IndsCoef,
integer, intent(in)  r 
)

Definition at line 525 of file BuildTensors.F90.

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 

◆ addtotensore()

recursive subroutine calctensore_list::addtotensore ( double complex, dimension(rts(r)), intent(in)  MomTenRec,
integer, dimension(4), intent(inout)  IndsCoef,
integer, intent(in)  r 
)

Definition at line 722 of file BuildTensors.F90.

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 

◆ addtotensorf()

recursive subroutine calctensorf_list::addtotensorf ( double complex, dimension(rts(r)), intent(in)  MomTenRec,
integer, dimension(5), intent(inout)  IndsCoef,
integer, intent(in)  r 
)

Definition at line 924 of file BuildTensors.F90.

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 

◆ addtotensorfuv()

recursive subroutine calctensorfuv_list::addtotensorfuv ( double complex, dimension(rts(r)), intent(in)  MomTenRec,
integer, dimension(5), intent(inout)  IndsCoef,
integer, intent(in)  r 
)

Definition at line 1095 of file BuildTensors.F90.

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 

◆ addtotensorg()

recursive subroutine calctensorg_list::addtotensorg ( double complex, dimension(rts(r)), intent(in)  MomTenRec,
integer, dimension(6), intent(inout)  IndsCoef,
integer  r 
)

Definition at line 1264 of file BuildTensors.F90.

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 

◆ addtotensortn()

recursive subroutine calctensortn_list::addtotensortn ( double complex, dimension(rts(r)), intent(in)  MomTenRec,
integer, intent(in)  ind,
integer, intent(in)  r 
)

Definition at line 1463 of file BuildTensors.F90.

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 

◆ addtotensortnuv()

recursive subroutine calctensortnuv_list::addtotensortnuv ( double complex, dimension(rts(r)), intent(in)  MomTenRec,
integer, intent(in)  ind,
integer, intent(in)  r 
)

Definition at line 1626 of file BuildTensors.F90.

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 
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