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

Go to the source code of this file.

Functions/Subroutines

subroutine attachline (pneigh, pnode, nterms, xil, etl, xn)
 

Function/Subroutine Documentation

◆ attachline()

subroutine attachline ( real*8, dimension(3,9), intent(in)  pneigh,
real*8, dimension(3), intent(in)  pnode,
integer, intent(in)  nterms,
real*8, intent(inout)  xil,
real*8, intent(inout)  etl,
real*8, dimension(3), intent(in)  xn 
)
20 !
21 ! returns the local coordinates in a face described by
22 ! the nodal coordinates in pneigh(1..3,*) of the intersection
23 ! with a straight line through point pnode(1..3) and
24 ! direction xn(1..3)
25 !
26  implicit none
27 !
28  integer nterms,i,j,k,imin,jmin,im,jm
29 !
30  real*8 pneigh(3,9),pnode(3),xi(-1:1,-1:1),
31  & et(-1:1,-1:1),distmin,d1,dist,xil,etl,xn(3)
32 !
33  intent(in) pneigh,nterms,xn,pnode
34 !
35  intent(inout) xil,etl
36 !
37  d1=1.d0
38 !
39  xi(0,0)=0.d0
40  et(0,0)=0.d0
41  call distattachline(xi(0,0),et(0,0),pneigh,pnode,dist,
42  & nterms,xn)
43  distmin=dist
44  imin=0
45  jmin=0
46 !
47  do k=1,6
48 !
49 ! initialisation
50 !
51  d1=d1/10.d0
52 !
53  do i=-1,1
54  do j=-1,1
55  if((i.eq.0).and.(j.eq.0)) cycle
56 !
57  xi(i,j)=xi(0,0)+i*d1
58  et(i,j)=et(0,0)+j*d1
59 !
60 ! check whether inside the (-1,1)x(-1,1) domain
61 !
62  if((xi(i,j).le.1.d0).and.
63  & (xi(i,j).ge.-1.d0).and.
64  & (et(i,j).le.1.d0).and.
65  & (et(i,j).ge.-1.d0)) then
66  call distattachline(xi(i,j),et(i,j),pneigh,pnode,
67  & dist,nterms,xn)
68 !
69 ! checking for smallest initial distance
70 !
71  if(dist.lt.distmin) then
72  distmin=dist
73  imin=i
74  jmin=j
75  endif
76  endif
77 !
78  enddo
79  enddo
80 !
81 ! minimizing the distance from the face to the node
82 !
83  do
84 !
85 ! exit if minimum found
86 !
87  if((imin.eq.0).and.(jmin.eq.0)) exit
88 !
89 ! new center of 3x3 matrix
90 !
91  xi(0,0)=xi(imin,jmin)
92  et(0,0)=et(imin,jmin)
93 !
94  im=imin
95  jm=jmin
96 !
97  imin=0
98  jmin=0
99 !
100  do i=-1,1
101  do j=-1,1
102  if((i+im.lt.-1).or.(i+im.gt.1).or.
103  & (j+jm.lt.-1).or.(j+jm.gt.1)) then
104 !
105  xi(i,j)=xi(0,0)+i*d1
106  et(i,j)=et(0,0)+j*d1
107 !
108 ! check whether inside the (-1,1)x(-1,1) domain
109 !
110  if((xi(i,j).le.1.d0).and.
111  & (xi(i,j).ge.-1.d0).and.
112  & (et(i,j).le.1.d0).and.
113  & (et(i,j).ge.-1.d0)) then
114  call distattachline(xi(i,j),et(i,j),pneigh,
115  & pnode,dist,nterms,xn)
116 !
117 ! check for new minimum
118 !
119  if(dist.lt.distmin) then
120  distmin=dist
121  imin=i
122  jmin=j
123  endif
124  endif
125 !
126  endif
127  enddo
128  enddo
129  enddo
130  enddo
131 !
132  if(nterms.eq.3) then
133  if(xi(0,0)+et(0,0).le.0.d0) then
134  xil=(xi(0,0)+1.d0)/2.d0
135  etl=(et(0,0)+1.d0)/2.d0
136  else
137  xil=(1.d0-et(0,0))/2.d0
138  etl=(1.d0-xi(0,0))/2.d0
139  endif
140  elseif(nterms.eq.4) then
141  xil=xi(0,0)
142  etl=et(0,0)
143  elseif(nterms.eq.6) then
144  if(xi(0,0)+et(0,0).le.0.d0) then
145  xil=(xi(0,0)+1.d0)/2.d0
146  etl=(et(0,0)+1.d0)/2.d0
147  else
148  xil=(1.d0-et(0,0))/2.d0
149  etl=(1.d0-xi(0,0))/2.d0
150  endif
151  elseif(nterms.eq.8) then
152  xil=xi(0,0)
153  etl=et(0,0)
154  endif
155 !
156  return
static double * dist
Definition: radflowload.c:42
subroutine distattachline(xig, etg, pneigh, pnode, dist, nterms, xn)
Definition: distattachline.f:21
Hosted by OpenAircraft.com, (Michigan UAV, LLC)