CalculiX  2.13
A Free Software Three-Dimensional Structural Finite Element Program
smalldist.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine smalldist (co, distmin, lakon, ipkon, kon, ne)
 

Function/Subroutine Documentation

◆ smalldist()

subroutine smalldist ( real*8, dimension(3,*)  co,
real*8  distmin,
character*8, dimension(*)  lakon,
integer, dimension(*)  ipkon,
integer, dimension(*)  kon,
integer  ne 
)
20 !
21  implicit none
22 !
23  character*8 lakon(*)
24 !
25  integer ipkon(*),kon(*),i,j,ne,n,j1,j2,indexe,neigh2(2,1),
26  & neigh4(2,6),neigh6(2,9),neigh8(2,12),neigh10(2,12),
27  & neigh15(2,18),neigh20(2,24)
28 !
29  real*8 dist,distmin,co(3,*)
30 !
31  neigh2=reshape((/1,2/),(/2,1/))
32  neigh4=reshape((/1,2,2,3,3,1,1,4,2,4,3,4/),(/2,6/))
33  neigh6=reshape((/1,2,2,3,3,1,4,5,5,6,6,4,1,4,2,5,3,6/),(/2,9/))
34  neigh8=reshape((/1,2,2,3,3,4,4,1,5,6,6,7,7,8,8,5,
35  & 1,5,2,6,3,7,4,8/),(/2,12/))
36  neigh10=reshape((/1,5,5,2,2,6,6,3,3,7,7,1,
37  & 1,8,8,4,2,9,9,4,3,10,10,4/),(/2,12/))
38  neigh15=reshape((/1,7,7,2,2,8,8,3,3,9,9,1,
39  & 4,10,10,5,5,11,11,6,6,12,12,4,
40  & 1,13,13,4,2,14,14,5,3,15,15,6/),(/2,18/))
41  neigh20=
42  & reshape((/1,9,9,2,2,10,10,3,3,11,11,4,4,12,12,1,
43  & 5,13,13,6,6,14,14,7,7,15,15,8,8,16,16,5,
44  & 1,17,17,5,2,18,18,6,3,19,19,7,4,20,20,8/),(/2,24/))
45 
46 ! determining the distance between nodes
47 !
48  distmin=1.d30
49 !
50  do i=1,ne
51  if(ipkon(i).lt.0) cycle
52  if(lakon(i)(1:1).eq.'F') cycle
53  indexe=ipkon(i)
54  if(lakon(i)(4:4).eq.'2') then
55  n=24
56  do j=1,n
57  j1=kon(indexe+neigh20(1,j))
58  j2=kon(indexe+neigh20(2,j))
59  dist=(co(1,j1)-co(1,j2))**2+
60  & (co(2,j1)-co(2,j2))**2+
61  & (co(3,j1)-co(3,j2))**2
62  distmin=min(dist,distmin)
63  enddo
64  elseif(lakon(i)(1:8).eq.'ESPRNGA1') then
65  n=1
66  do j=1,n
67  j1=kon(indexe+neigh2(1,j))
68  j2=kon(indexe+neigh2(2,j))
69  dist=(co(1,j1)-co(1,j2))**2+
70  & (co(2,j1)-co(2,j2))**2+
71  & (co(3,j1)-co(3,j2))**2
72  distmin=min(dist,distmin)
73  enddo
74  elseif(lakon(i)(4:4).eq.'8') then
75  n=12
76  do j=1,n
77  j1=kon(indexe+neigh8(1,j))
78  j2=kon(indexe+neigh8(2,j))
79  dist=(co(1,j1)-co(1,j2))**2+
80  & (co(2,j1)-co(2,j2))**2+
81  & (co(3,j1)-co(3,j2))**2
82  distmin=min(dist,distmin)
83  enddo
84  elseif(lakon(i)(4:4).eq.'4') then
85  n=6
86  do j=1,n
87  j1=kon(indexe+neigh4(1,j))
88  j2=kon(indexe+neigh4(2,j))
89  dist=(co(1,j1)-co(1,j2))**2+
90  & (co(2,j1)-co(2,j2))**2+
91  & (co(3,j1)-co(3,j2))**2
92  distmin=min(dist,distmin)
93  enddo
94  elseif(lakon(i)(4:5).eq.'10') then
95  n=12
96  do j=1,n
97  j1=kon(indexe+neigh10(1,j))
98  j2=kon(indexe+neigh10(2,j))
99  dist=(co(1,j1)-co(1,j2))**2+
100  & (co(2,j1)-co(2,j2))**2+
101  & (co(3,j1)-co(3,j2))**2
102  distmin=min(dist,distmin)
103  enddo
104  elseif(lakon(i)(4:4).eq.'6') then
105  n=9
106  do j=1,n
107  j1=kon(indexe+neigh6(1,j))
108  j2=kon(indexe+neigh6(2,j))
109  dist=(co(1,j1)-co(1,j2))**2+
110  & (co(2,j1)-co(2,j2))**2+
111  & (co(3,j1)-co(3,j2))**2
112  distmin=min(dist,distmin)
113  enddo
114  elseif(lakon(i)(4:5).eq.'15') then
115  n=18
116  do j=1,n
117  j1=kon(indexe+neigh15(1,j))
118  j2=kon(indexe+neigh15(2,j))
119  dist=(co(1,j1)-co(1,j2))**2+
120  & (co(2,j1)-co(2,j2))**2+
121  & (co(3,j1)-co(3,j2))**2
122  distmin=min(dist,distmin)
123  enddo
124  else
125  cycle
126  endif
127  enddo
128 !
129 c distmin=dsqrt(distmin)*1.0e-03
130  distmin=dsqrt(distmin)*1.0e-04
131 !
132  return
#define min(a, b)
Definition: cascade.c:31
static double * dist
Definition: radflowload.c:42
Hosted by OpenAircraft.com, (Michigan UAV, LLC)