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

Go to the source code of this file.

Functions/Subroutines

subroutine attach (pneigh, pnode, nterms, ratio, dist, xil, etl)
 

Function/Subroutine Documentation

◆ attach()

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