xref: /phasta/phSolver/compressible/d2wall.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen
2*59599516SKenneth E. Jansen        subroutine D2wall (dwall,       x)
3*59599516SKenneth E. Jansenc
4*59599516SKenneth E. Jansenc----------------------------------------------------------------------
5*59599516SKenneth E. Jansenc
6*59599516SKenneth E. Jansenc This routine calculates the distance to the nearest wall for all
7*59599516SKenneth E. Jansenc field points.
8*59599516SKenneth E. Jansenc
9*59599516SKenneth E. Jansenc
10*59599516SKenneth E. Jansenc Zdenek Johan, Winter 1991.  (Fortran 90)
11*59599516SKenneth E. Jansenc----------------------------------------------------------------------
12*59599516SKenneth E. Jansenc
13*59599516SKenneth E. Jansen        use pointer_data
14*59599516SKenneth E. Jansen        include "common.h"
15*59599516SKenneth E. Jansen        include "mpif.h"
16*59599516SKenneth E. Jansenc
17*59599516SKenneth E. Jansen        dimension x(numnp,nsd),        dwall(numnp),
18*59599516SKenneth E. Jansen     &            nwall(numpe),        idisp(numpe)
19*59599516SKenneth E. Jansen
20*59599516SKenneth E. Jansenc
21*59599516SKenneth E. Jansenc  Allocatable arrays
22*59599516SKenneth E. Jansenc
23*59599516SKenneth E. Jansen        allocatable xwi(:,:,:),  xw(:,:,:)
24*59599516SKenneth E. Jansen
25*59599516SKenneth E. Jansenc
26*59599516SKenneth E. Jansenc.... Count the welts (wall-elements)
27*59599516SKenneth E. Jansenc
28*59599516SKenneth E. Jansen        nwalli=0
29*59599516SKenneth E. Jansen        do iblk = 1, nelblb ! loop over boundary elt blocks
30*59599516SKenneth E. Jansen           npro = lcblkb(1,iblk+1) - lcblkb(1,iblk)
31*59599516SKenneth E. Jansen           do j = 1, npro
32*59599516SKenneth E. Jansen              if(miBCB(iblk)%p(j,5).eq.2) nwalli=nwalli+1
33*59599516SKenneth E. Jansen           enddo
34*59599516SKenneth E. Jansen        enddo
35*59599516SKenneth E. Jansenc
36*59599516SKenneth E. Jansenc.... Create wallnode-coord list for welts on processor
37*59599516SKenneth E. Jansenc
38*59599516SKenneth E. Jansen        if (nwalli.eq.0) nwalli=1  !  patch for mpi's lack of imagination
39*59599516SKenneth E. Jansen        allocate (xwi(nwalli,nenb+1,nsd))
40*59599516SKenneth E. Jansen        xwi = 1.0d32
41*59599516SKenneth E. Jansen        xwi(:,nenb+1,:)=zero
42*59599516SKenneth E. Jansen        nwalli = 0
43*59599516SKenneth E. Jansen        do iblk = 1, nelblb ! loop over boundary elt blocks
44*59599516SKenneth E. Jansenc
45*59599516SKenneth E. Jansen           iel    = lcblkb(1,iblk)
46*59599516SKenneth E. Jansen           nenbl  = lcblkb(6,iblk) ! no. of vertices per bdry. face
47*59599516SKenneth E. Jansen           npro   = lcblkb(1,iblk+1) - iel
48*59599516SKenneth E. Jansenc
49*59599516SKenneth E. Jansen           do j = 1, npro ! loop over belts in this blk
50*59599516SKenneth E. Jansen              if(miBCB(iblk)%p(j,5).eq.2) then
51*59599516SKenneth E. Jansen                 nwalli=nwalli+1
52*59599516SKenneth E. Jansenc assemble local coord list
53*59599516SKenneth E. Jansen                 do node = 1, nenbl
54*59599516SKenneth E. Jansen                    xwi(nwalli,node,1:3)=x(mienb(iblk)%p(j,node),:)
55*59599516SKenneth E. Jansen                 enddo
56*59599516SKenneth E. Jansen                 do node = 1, nenbl
57*59599516SKenneth E. Jansen                    xwi(nwalli,nenb+1,:)=xwi(nwalli,nenb+1,:)
58*59599516SKenneth E. Jansen     &                                   +xwi(nwalli,node,:)
59*59599516SKenneth E. Jansen                 enddo
60*59599516SKenneth E. Jansen                 xwi(nwalli,nenb+1,:)=xwi(nwalli,nenb+1,:)/nenbl
61*59599516SKenneth E. Jansenc
62*59599516SKenneth E. Jansen              endif
63*59599516SKenneth E. Jansen           enddo          ! loop over belts in this blk
64*59599516SKenneth E. Jansenc
65*59599516SKenneth E. Jansen        enddo               ! loop over boundary elt blocks
66*59599516SKenneth E. Jansenc
67*59599516SKenneth E. Jansen        if (nwalli.eq.0) xwi=1.0e32 ! fix for mpi's lack of imagination
68*59599516SKenneth E. Jansen        if (nwalli.eq.0) nwalli=1  !  patch for mpi's lack of imagination
69*59599516SKenneth E. Jansenc
70*59599516SKenneth E. Jansenc  Pool "number of welts" info from all processors
71*59599516SKenneth E. Jansenc
72*59599516SKenneth E. Jansen        if (numpe.gt.1) then
73*59599516SKenneth E. Jansen           call MPI_ALLGATHER(nwalli,1,MPI_INTEGER,nwall,1,
74*59599516SKenneth E. Jansen     &          MPI_INTEGER,MPI_COMM_WORLD,ierr)
75*59599516SKenneth E. Jansen           nwallt=0
76*59599516SKenneth E. Jansen           do j=1,numpe
77*59599516SKenneth E. Jansen              nwallt=nwallt+nwall(j)
78*59599516SKenneth E. Jansen           enddo
79*59599516SKenneth E. Jansen        else
80*59599516SKenneth E. Jansen           nwall=nwalli
81*59599516SKenneth E. Jansen           nwallt=nwalli
82*59599516SKenneth E. Jansen        endif
83*59599516SKenneth E. Jansenc
84*59599516SKenneth E. Jansenc  Make all-processor wallnode-coord collage
85*59599516SKenneth E. Jansenc
86*59599516SKenneth E. Jansen        allocate (xw(nwallt,nenb+1,nsd))
87*59599516SKenneth E. Jansen        if (numpe.gt.1) then
88*59599516SKenneth E. Jansen           idisp(:)=0
89*59599516SKenneth E. Jansen           do j=2,numpe
90*59599516SKenneth E. Jansen              idisp(j)=idisp(j-1)+nwall(j-1)
91*59599516SKenneth E. Jansen           enddo
92*59599516SKenneth E. Jansen           do j=1,nenb+1
93*59599516SKenneth E. Jansen              do k=1,nsd
94*59599516SKenneth E. Jansen                 call MPI_ALLGATHERV(xwi(:,j,k),nwalli,
95*59599516SKenneth E. Jansen     &                MPI_DOUBLE_PRECISION,xw(:,j,k),nwall,idisp,
96*59599516SKenneth E. Jansen     &                MPI_DOUBLE_PRECISION,MPI_COMM_WORLD,ierr)
97*59599516SKenneth E. Jansen              enddo
98*59599516SKenneth E. Jansen           enddo
99*59599516SKenneth E. Jansen        else
100*59599516SKenneth E. Jansen           xw=xwi
101*59599516SKenneth E. Jansen        endif
102*59599516SKenneth E. Jansen
103*59599516SKenneth E. Jansenc
104*59599516SKenneth E. Jansenc  For each point, calculate distance to centroid;
105*59599516SKenneth E. Jansenc  save the distance in this node's position of dwall if it's
106*59599516SKenneth E. Jansenc  shorter than currently stored distance
107*59599516SKenneth E. Jansenc
108*59599516SKenneth E. Jansen        dwall=1.0e32
109*59599516SKenneth E. Jansen        do i=1,numnp
110*59599516SKenneth E. Jansen           do j=1, nwallt
111*59599516SKenneth E. Jansen              do k=1,nenb+1
112*59599516SKenneth E. Jansen                 distance =  ( x(i,1) - xw(j,k,1) )**2
113*59599516SKenneth E. Jansen     &                      +( x(i,2) - xw(j,k,2) )**2
114*59599516SKenneth E. Jansen     &                      +( x(i,3) - xw(j,k,3) )**2
115*59599516SKenneth E. Jansen                 if ( dwall(i).gt.distance ) dwall(i) = distance
116*59599516SKenneth E. Jansen              enddo
117*59599516SKenneth E. Jansen           enddo
118*59599516SKenneth E. Jansen        enddo
119*59599516SKenneth E. Jansen        dwall=sqrt(dwall)
120*59599516SKenneth E. Jansenc
121*59599516SKenneth E. Jansen        deallocate(xwi)
122*59599516SKenneth E. Jansen        deallocate(xw)
123*59599516SKenneth E. Jansenc
124*59599516SKenneth E. Jansen        return
125*59599516SKenneth E. Jansen        end
126*59599516SKenneth E. Jansen
127*59599516SKenneth E. Jansen
128*59599516SKenneth E. Jansen
129