xref: /phasta/phSolver/common/sonfath.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen      subroutine sonfath (ncorp, vnumnpl, ifath, nsons, tfath,
2*59599516SKenneth E. Jansen     &                    numpe, numnp,   maxnode)
3*59599516SKenneth E. Jansen
4*59599516SKenneth E. Jansenc     This program writes to geom.dat.<1..numpe> the number sons and
5*59599516SKenneth E. Jansenc     fathers required for averaging in the z-direction (Turbulent
6*59599516SKenneth E. Jansenc     Boundary Layer), (x,y)-space (Turbulent Channel Flow),
7*59599516SKenneth E. Jansenc     or (x,y,z)-space (Isotropic Turbulence)
8*59599516SKenneth E. Jansen
9*59599516SKenneth E. Jansen      integer inform,n1m,n2m,n3m,numpe,numpe1,numnp,maxnode, tfath,
10*59599516SKenneth E. Jansen     &     counter,nc8,nc4,node,ifil,itmp
11*59599516SKenneth E. Jansen      integer gfath(0:numnp),ifath(numpe,maxnode),ncorp(numpe,maxnode),
12*59599516SKenneth E. Jansen     &     nsons(tfath),dummy(numnp),iper(numnp), vnumnpl(numpe)
13*59599516SKenneth E. Jansen      integer, allocatable :: imap(:),invmap(:)
14*59599516SKenneth E. Jansen      allocate (imap(numnp))
15*59599516SKenneth E. Jansen      allocate (invmap(numnp))
16*59599516SKenneth E. Jansen
17*59599516SKenneth E. Jansen      write(*,*)'Enter number of nodes in each direction'
18*59599516SKenneth E. Jansen      read(*,*)n1m,n2m,n3m
19*59599516SKenneth E. Jansen
20*59599516SKenneth E. Jansen      write(*,*)'For Z averaging enter 1 For X-Y averaging enter 2'
21*59599516SKenneth E. Jansen      write(*,*)'For X-Y-Z averaging enter 3  W=0 SPEBC enter 4'
22*59599516SKenneth E. Jansen      write(*,*)'fully local model 5'
23*59599516SKenneth E. Jansen      read(*,*)inform
24*59599516SKenneth E. Jansen
25*59599516SKenneth E. Jansenc
26*59599516SKenneth E. Jansenc.... read and generate map
27*59599516SKenneth E. Jansenc
28*59599516SKenneth E. Jansen      kmap = 15
29*59599516SKenneth E. Jansen
30*59599516SKenneth E. Jansen      open(kmap,file='geom.kmap',status='old',err=183)
31*59599516SKenneth E. Jansen
32*59599516SKenneth E. Jansen      do i=1,numnp
33*59599516SKenneth E. Jansen         read(kmap,*,err=183) iorig
34*59599516SKenneth E. Jansen         imap(i)=iorig+1  ! imap takes arguement of original and returns new
35*59599516SKenneth E. Jansen         invmap(iorig+1)=i ! invmap takes arguemnt of new and returns original
36*59599516SKenneth E. Jansen      enddo
37*59599516SKenneth E. Jansen      goto 184
38*59599516SKenneth E. Jansen 183  write(*,*) 'geom.kmap not read properly...assuming inform=5'
39*59599516SKenneth E. Jansen      inform=5
40*59599516SKenneth E. Jansen 184  continue
41*59599516SKenneth E. Jansen
42*59599516SKenneth E. Jansen      nsoncor=1
43*59599516SKenneth E. Jansen      if(inform.eq.4) then
44*59599516SKenneth E. Jansen         inform=1
45*59599516SKenneth E. Jansen         nsoncor=0
46*59599516SKenneth E. Jansen      endif
47*59599516SKenneth E. Jansen
48*59599516SKenneth E. Jansen      gfath(:) = 0
49*59599516SKenneth E. Jansen
50*59599516SKenneth E. Jansen      nc4=n3m*(n2m-1)+1
51*59599516SKenneth E. Jansen      nc8=n3m*n2m*n1m+1-n3m
52*59599516SKenneth E. Jansen
53*59599516SKenneth E. Jansen      if (numpe.eq.1)then
54*59599516SKenneth E. Jansen         numnp   = n1m*n2m*n3m
55*59599516SKenneth E. Jansen         maxnode = numnp
56*59599516SKenneth E. Jansen      endif
57*59599516SKenneth E. Jansen
58*59599516SKenneth E. Jansen      if (inform.eq.1) then ! For averaging over Z - lines
59*59599516SKenneth E. Jansen
60*59599516SKenneth E. Jansen         do i = 1, tfath
61*59599516SKenneth E. Jansen            nsons(i) = n3m-nsoncor
62*59599516SKenneth E. Jansen         enddo
63*59599516SKenneth E. Jansen
64*59599516SKenneth E. Jansen         counter = 1
65*59599516SKenneth E. Jansen         do i = 1, nc8, n3m
66*59599516SKenneth E. Jansen            do j = i, i+n3m-1
67*59599516SKenneth E. Jansen               gfath(j) = counter
68*59599516SKenneth E. Jansen            enddo
69*59599516SKenneth E. Jansen            counter = counter + 1
70*59599516SKenneth E. Jansen         enddo
71*59599516SKenneth E. Jansen
72*59599516SKenneth E. Jansen         if(numpe.gt.1)then
73*59599516SKenneth E. Jansen            do j = 1, numpe
74*59599516SKenneth E. Jansen               do i = 1, maxnode
75*59599516SKenneth E. Jansen                  ifath(j,i) = gfath(imap(ncorp(j,i)))
76*59599516SKenneth E. Jansen               enddo
77*59599516SKenneth E. Jansen            enddo
78*59599516SKenneth E. Jansen         else
79*59599516SKenneth E. Jansen            do i = 1, maxnode
80*59599516SKenneth E. Jansen               ifath(1,i) = gfath(imap(i))
81*59599516SKenneth E. Jansen            enddo
82*59599516SKenneth E. Jansen         endif
83*59599516SKenneth E. Jansen
84*59599516SKenneth E. Jansen      endif
85*59599516SKenneth E. Jansen
86*59599516SKenneth E. Jansen      if (inform.eq.2) then ! For averaging over X-Z planes
87*59599516SKenneth E. Jansen
88*59599516SKenneth E. Jansen         do i = 1, tfath
89*59599516SKenneth E. Jansen            nsons(i) = n1m*n3m - (n1m+n3m-1)
90*59599516SKenneth E. Jansen         enddo
91*59599516SKenneth E. Jansen
92*59599516SKenneth E. Jansen         counter = 1
93*59599516SKenneth E. Jansen         do j = 1, nc4, n3m
94*59599516SKenneth E. Jansen            do i = 1, n1m
95*59599516SKenneth E. Jansen               do k = j, j+n3m-1
96*59599516SKenneth E. Jansen                  node = (i-1)*n2m*n3m + k
97*59599516SKenneth E. Jansen                  gfath(node) = counter
98*59599516SKenneth E. Jansen               enddo
99*59599516SKenneth E. Jansen            enddo
100*59599516SKenneth E. Jansen            counter = counter + 1
101*59599516SKenneth E. Jansen         enddo
102*59599516SKenneth E. Jansen
103*59599516SKenneth E. Jansen         if(numpe.gt.1)then
104*59599516SKenneth E. Jansen            do j = 1, numpe
105*59599516SKenneth E. Jansen               do i = 1, maxnode
106*59599516SKenneth E. Jansen                  ifath(j,i) = gfath(imap(ncorp(j,i)))
107*59599516SKenneth E. Jansen               enddo
108*59599516SKenneth E. Jansen            enddo
109*59599516SKenneth E. Jansen         else
110*59599516SKenneth E. Jansen            do i = 1, maxnode
111*59599516SKenneth E. Jansen               ifath(1,i) = gfath(imap(i))
112*59599516SKenneth E. Jansen            enddo
113*59599516SKenneth E. Jansen         endif
114*59599516SKenneth E. Jansen
115*59599516SKenneth E. Jansen      endif
116*59599516SKenneth E. Jansen
117*59599516SKenneth E. Jansen      if( inform.eq.3 )then !For averaging over X-Y-Z
118*59599516SKenneth E. Jansen
119*59599516SKenneth E. Jansen         open (unit=96,FILE='geom.iper',form='unformatted',
120*59599516SKenneth E. Jansen     &        status='old')
121*59599516SKenneth E. Jansen         read(96) (iper(l),l=1,numnp)
122*59599516SKenneth E. Jansen         close(96)
123*59599516SKenneth E. Jansen
124*59599516SKenneth E. Jansen         dummy(:) = 1
125*59599516SKenneth E. Jansen
126*59599516SKenneth E. Jansen         do i = 1, numnp
127*59599516SKenneth E. Jansen            if(iper(i).eq.0)then
128*59599516SKenneth E. Jansen               iper(i)=i
129*59599516SKenneth E. Jansen            endif
130*59599516SKenneth E. Jansen         enddo
131*59599516SKenneth E. Jansen         do j = 1, numnp
132*59599516SKenneth E. Jansen            i = iper(j)
133*59599516SKenneth E. Jansen            if (i .ne. j) then
134*59599516SKenneth E. Jansen               dummy(j) = 0
135*59599516SKenneth E. Jansen            endif
136*59599516SKenneth E. Jansen         enddo
137*59599516SKenneth E. Jansen         itmp=sum(dummy)
138*59599516SKenneth E. Jansen         nsons(tfath) = itmp
139*59599516SKenneth E. Jansen         write(*,*)'nsons(1)=',nsons(1)
140*59599516SKenneth E. Jansen         ifath=1
141*59599516SKenneth E. Jansen      endif
142*59599516SKenneth E. Jansen
143*59599516SKenneth E. Jansen      if(inform.eq.5) then      ! local clipping
144*59599516SKenneth E. Jansen
145*59599516SKenneth E. Jansen         do i = 1, numpe
146*59599516SKenneth E. Jansen            ifil = 99
147*59599516SKenneth E. Jansen
148*59599516SKenneth E. Jansen            nsons=1
149*59599516SKenneth E. Jansen            do j=1,vnumnpl(i)
150*59599516SKenneth E. Jansen               ifath(i,j)=j
151*59599516SKenneth E. Jansen            enddo
152*59599516SKenneth E. Jansen
153*59599516SKenneth E. Jansen         enddo
154*59599516SKenneth E. Jansen      endif
155*59599516SKenneth E. Jansen
156*59599516SKenneth E. Jansen      return
157*59599516SKenneth E. Jansen      end
158