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

Go to the source code of this file.

Functions/Subroutines

subroutine distattachline (xig, etg, pneigh, pnode, dist, nterms, xn)
 

Function/Subroutine Documentation

◆ distattachline()

subroutine distattachline ( real*8  xig,
real*8  etg,
real*8, dimension(3,*)  pneigh,
real*8, dimension(3)  pnode,
real*8, intent(inout)  dist,
integer  nterms,
real*8, dimension(3)  xn 
)
21 !
22 ! calculates the distance between a straight line through the node
23 ! with coordinates in "pnode" and direction vector "xn" and
24 ! the node with local coordinates xig and etg
25 ! in a face described by "nterms" nodes with coordinates
26 ! in "pneigh"
27 !
28  implicit none
29 !
30  integer nterms,i,j
31 !
32  real*8 ratio(8),pneigh(3,*),pnode(3),dist,xi,et,xig,etg,p(3),
33  & xn(3),coeff,etm2,xim2,etm,xim,etp,xip,a2,et2,xi2,a
34 !
35  intent(in)xig,etg,pneigh,pnode,nterms,xn
36 !
37  intent(inout) dist
38 !
39  if(nterms.eq.3) then
40  if(xig+etg.le.0.d0) then
41  ratio(2)=(xig+1.d0)/2.d0
42  ratio(3)=(etg+1.d0)/2.d0
43  else
44  ratio(2)=(1.d0-etg)/2.d0
45  ratio(3)=(1.d0-xig)/2.d0
46  endif
47  ratio(1)=1.d0-ratio(2)-ratio(3)
48  elseif(nterms.eq.4) then
49  xip=(1.d0+xig)/4.d0
50  xim=(1.d0-xig)/4.d0
51  etp=1.d0+etg
52  etm=1.d0-etg
53  ratio(1)=xim*etm
54  ratio(2)=xip*etm
55  ratio(3)=xip*etp
56  ratio(4)=xim*etp
57  elseif(nterms.eq.6) then
58  if(xig+etg.le.0.d0) then
59  xi=(xig+1.d0)/2.d0
60  et=(etg+1.d0)/2.d0
61  else
62  xi=(1.d0-etg)/2.d0
63  et=(1.d0-xig)/2.d0
64  endif
65  a=1.d0-xi-et
66  a2=2.d0*a
67  xi2=2.d0*xi
68  et2=2.d0*et
69  ratio(1)=a*(a2-1.d0)
70  ratio(2)=xi*(xi2-1.d0)
71  ratio(3)=et*(et2-1.d0)
72  ratio(4)=xi2*a2
73  ratio(5)=xi2*et2
74  ratio(6)=et2*a2
75  elseif(nterms.eq.8) then
76  xip=1.d0+xig
77  xim=1.d0-xig
78  xim2=xip*xim/2.d0
79  etp=1.d0+etg
80  etm=1.d0-etg
81  etm2=etp*etm/2.d0
82  ratio(5)=xim2*etm
83  ratio(6)=xip*etm2
84  ratio(7)=xim2*etp
85  ratio(8)=xim*etm2
86  xim=xim/4.d0
87  xip=xip/4.d0
88  ratio(1)=xim*etm*(-xig-etp)
89  ratio(2)=xip*etm*(xig-etp)
90  ratio(3)=xip*etp*(xig-etm)
91  ratio(4)=xim*etp*(-xig-etm)
92  else
93  write(*,*) '*ERROR in distattach: case with ',nterms
94  write(*,*) ' terms is not covered'
95  call exit(201)
96  endif
97 !
98 ! calculating the position in the face
99 !
100  do i=1,3
101  p(i)=0.d0
102  do j=1,nterms
103  p(i)=p(i)+ratio(j)*pneigh(i,j)
104  enddo
105  enddo
106 !
107 ! calculating the distance
108 !
109  coeff=0.0
110  do i=1,3
111  coeff=coeff+xn(i)*(p(i)-pnode(i))
112  enddo
113  dist=(p(1)-pnode(1)-coeff*xn(1))**2+(p(2)-pnode(2)-
114  & coeff*xn(2))**2+(p(3)-pnode(3)-coeff*xn(3))**2
115 !
116  return
static double * dist
Definition: radflowload.c:42
Hosted by OpenAircraft.com, (Michigan UAV, LLC)