summaryrefslogtreecommitdiff
path: root/Dragon/src/g2s_generatingTrack.f90
blob: 14501ee0d5a935b7ebae4496e572e7e1b68a8a0e (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
!
!-----------------------------------------------------------------------
!
!Purpose:
! Generate data relative to gigognes originating from nodes and generate
! tracking indices assigned to them.
!
!Copyright:
! Copyright (C) 2001 Ecole Polytechnique de Montreal.
! This library is free software; you can redistribute it and/or
! modify it under the terms of the GNU Lesser General Public
! License as published by the Free Software Foundation; either
! version 2.1 of the License, or (at your option) any later version.
!
!Author(s):
! G. Civario (CS-SI)
!
!Comments:
! Ce module possede trois fonctions:
!  - generateTrack : fonction d'entree du module
!  - calculFinalMerge : calcul de la numerotation type dragon
!  - ltVec : fonction d'ordre specifique pour des vecteurs d'entiers de tailles
!            differentes
!
module track

    use cellulePlaced
    use segArc
    use ptNodes
    implicit none

contains
  subroutine generateTrack(szP,szSa,nbNode,lgMaxGig,gig,merg)
    integer,intent(in) :: szP,szSa,nbNode,lgMaxGig
    integer,dimension(lgMaxGig*nbNode),intent(out) :: gig
    integer,dimension(nbNode),intent(inout) :: merg
    
    integer :: i,lgMaxMrg,s,d,lgMaxGigTest
    integer,dimension(:,:),allocatable :: mrgMat
    
    lgMaxGigTest = 0
    lgMaxMrg = 0
    do i = 1,szP
       lgMaxGigTest = max(lgMaxGigTest,size(tabCellulePlaced(i)%gig))
       lgMaxMrg = max(lgMaxMrg,size(tabCellulePlaced(i)%mrg))
    end do
    if(lgMaxGigTest /= lgMaxGig) call XABORT('g2s_generatingTrack: lgMax error')
    lgMaxMrg = lgMaxMrg + 1
    allocate(mrgMat(lgMaxMrg,nbNode),stat=alloc_ok)
    if (alloc_ok /= 0) call XABORT("G2S: generateTrack => allocation pb")
    gig = 0
    mrgMat = 0
    do i = 1,szSa
       if (tabSegArc(i)%nodeg>0.and.tabSegArc(i)%indCellPg>0) then
          s = size(tabCellulePlaced(tabSegArc(i)%indCellPg)%gig)
          d = (tabSegArc(i)%nodeg-1)*lgMaxGig+1
          gig(d:d+s-1) = tabCellulePlaced(tabSegArc(i)%indCellPg)%gig(1:s)
          s = size(tabCellulePlaced(tabSegArc(i)%indCellPg)%mrg)
          d = tabSegArc(i)%nodeg
          mrgMat(1:s,d) = tabCellulePlaced(tabSegArc(i)%indCellPg)%mrg(1:s)
          mrgMat(s+1,d) = merg(d)
       end if
       if (tabSegArc(i)%noded>0.and.tabSegArc(i)%indCellPd>0) then
          s = size(tabCellulePlaced(tabSegArc(i)%indCellPd)%gig)
          d = (tabSegArc(i)%noded-1)*lgMaxGig+1
          gig(d:d+s-1) = tabCellulePlaced(tabSegArc(i)%indCellPd)%gig(1:s)
          s = size(tabCellulePlaced(tabSegArc(i)%indCellPd)%mrg)
          d = tabSegArc(i)%noded
          mrgMat(1:s,d) = tabCellulePlaced(tabSegArc(i)%indCellPd)%mrg(1:s)
          mrgMat(s+1,d) = merg(d)
       end if
    end do
    call calculFinalMerge(mrgMat,merg)
    deallocate(mrgMat)
  end subroutine generateTrack

  subroutine calculFinalMerge(inMat,outVec)
    integer,dimension(:,:),intent(in) :: inMat
    integer,dimension(:),intent(out)  :: outVec
    
    integer :: i,j,d1,d2,maxD2
    logical :: found,sorted
    integer,dimension(:),allocatable   :: tmpVec
    integer,dimension(:,:),allocatable :: workMat
    
    d1 = size(inMat,1) !profondeur max de gigogne + 1
    d2 = size(inMat,2) !nombre de nodes
    maxD2 = 0
    allocate(workMat(d1,d2))
    allocate(tmpVec(d1),stat=alloc_ok)
    if (alloc_ok /= 0) call XABORT("G2S: calculFinalMerge => allocation pb")
    workMat(:d1,:d2) = 0
    !remplissage de workMat a l'aide d'occurences uniques de lignes de inMat
    do i = 1,d2
       found = .false.
       do j = 1,maxD2
          if (all(workMat(:d1,j)==inMat(:d1,i))) then
             found = .true.
             exit
          end if
       end do
       if (.not. found) then
          maxD2 = maxD2 + 1
          workMat(:d1,maxD2) = inMat(:d1,i)
       end if
    end do
    !classement des lignes de workMat en ordre croissant par bubble-sort
    do j = maxD2,2,-1
       sorted = .true.
       do i = 1,j-1
          if (ltVec(workMat(:d1,i+1),workMat(:d1,i))) then
             tmpVec(:d1)      = workMat(:d1,i+1)
             workMat(:d1,i+1) = workMat(:d1,i)
             workMat(:d1,i)   = tmpVec(:d1)
             sorted = .false.
          end if
       end do
       if (sorted) exit
    end do
    !remplissage de outVec en fonction de l'egalite entre les lignes de workMat
    !et inMat, apres le classement
    do i = 1,d2
       do j = 1,maxD2
          if (all(workMat(:d1,j)==inMat(:d1,i))) then
             outVec(i) = j
             exit
          end if
       end do
    end do
    
    deallocate(workMat,tmpVec)
  end subroutine calculFinalMerge
  
  function ltVec(v1,v2)
    integer,dimension(:),intent(in) :: v1,v2
    logical :: ltVec
    integer :: i
    
    do i = 1,size(v1)
       if (v1(i) < v2(i)) then
          ltVec = .true.
          return
       else if (v1(i) > v2(i)) then
          ltVec = .false.
          return
       end if
    end do
    ltVec = .false.
  end function ltVec
end module track