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

Go to the source code of this file.

Functions/Subroutines

subroutine attach_3d (pneigh, pnode, nterms, ratio, dist, xil, etl, zel)
 

Function/Subroutine Documentation

◆ attach_3d()

subroutine attach_3d ( real*8, dimension(3,20), intent(in)  pneigh,
real*8, dimension(3), intent(inout)  pnode,
integer, intent(in)  nterms,
real*8, dimension(20), intent(inout)  ratio,
real*8, intent(inout)  dist,
real*8, intent(inout)  xil,
real*8, intent(inout)  etl,
real*8, intent(inout)  zel 
)
20 !
21 ! ataches node with coordinates in "pnode" to the face containing
22 ! "nterms" nodes with coordinates in field "pneigh" (nterms < 21).
23 ! cave: the coordinates are stored in pneigh(1..3,*)
24 !
25  implicit none
26 !
27  integer nterms,i,j,k,imin,jmin,kmin,kflag,n,iy
28 !
29  real*8 ratio(20),pneigh(3,20),pnode(3),dummy,
30  & a(-1:1,-1:1,-1:1),xi(-1:1,-1:1,-1:1),et(-1:1,-1:1,-1:1),p(3),
31  & aold(-1:1,-1:1,-1:1),ze(-1:1,-1:1,-1:1),zeold(-1:1,-1:1,-1:1),
32  & xiold(-1:1,-1:1,-1:1),etold(-1:1,-1:1,-1:1),distmin,xiopt,etopt,
33  & d1,d2,d3,d4,dist,xil,etl,zel,zeopt,dx(3),al
34 !
35  intent(in) pneigh,nterms
36 !
37  intent(inout) xil,etl,zel,dist,pnode,ratio
38 !
39  kflag=1
40  n=3
41 !
42 c d1=0.25d0
43 c d2=3.125d-2
44 c d3=3.9063d-3
45 c d4=1.d-3
46  d1=1.d-2
47  d2=1.d-4
48  d3=1.d-6
49  d4=1.d-7
50 !
51 ! initialisation
52 !
53  do i=-1,1
54  do j=-1,1
55  do k=-1,1
56  xi(i,j,k)=i*d1
57  et(i,j,k)=j*d1
58  ze(i,j,k)=k*d1
59  call distattach_3d(xi(i,j,k),et(i,j,k),ze(i,j,k),pneigh,
60  & pnode,a(i,j,k),p,ratio,nterms)
61 c write(*,'(3(1x,i5),4(1x,e11.4))')
62 c & i,j,k,xi(i,j,k),et(i,j,k),ze(i,j,k),a(i,j,k)
63  enddo
64  enddo
65  enddo
66 !
67 ! minimizing the distance from the face to the node
68 !
69  do
70  distmin=a(0,0,0)
71  imin=0
72  jmin=0
73  kmin=0
74 c write(*,*) '1 ',xi(0,0,0),et(0,0,0),
75 c & ze(0,0,0),distmin
76  do i=-1,1
77  do j=-1,1
78  do k=-1,1
79  if(a(i,j,k).lt.distmin) then
80  distmin=a(i,j,k)
81  imin=i
82  jmin=j
83  kmin=k
84  endif
85  enddo
86  enddo
87  enddo
88 c write(*,*) '1 ',xi(imin,jmin,kmin),et(imin,jmin,kmin),
89 c & ze(imin,jmin,kmin),distmin
90 cc write(*,*) '1 ',xi(0,0,0),et(0,0,0),ze(0,0,0)
91 !
92 ! exit if minimum found
93 !
94  if((imin.eq.0).and.(jmin.eq.0).and.(kmin.eq.0)) exit
95 !
96  do i=-1,1
97  do j=-1,1
98  do k=-1,1
99  aold(i,j,k)=a(i,j,k)
100  xiold(i,j,k)=xi(i,j,k)
101  etold(i,j,k)=et(i,j,k)
102  zeold(i,j,k)=ze(i,j,k)
103  enddo
104  enddo
105  enddo
106 !
107  do i=-1,1
108  do j=-1,1
109  do k=-1,1
110  if((i+imin.ge.-1).and.(i+imin.le.1).and.
111  & (j+jmin.ge.-1).and.(j+jmin.le.1).and.
112  & (k+kmin.ge.-1).and.(k+kmin.le.1)) then
113  a(i,j,k)=aold(i+imin,j+jmin,k+kmin)
114  xi(i,j,k)=xiold(i+imin,j+jmin,k+kmin)
115  et(i,j,k)=etold(i+imin,j+jmin,k+kmin)
116  ze(i,j,k)=zeold(i+imin,j+jmin,k+kmin)
117  else
118  xi(i,j,k)=xi(i,j,k)+imin*d1
119  et(i,j,k)=et(i,j,k)+jmin*d1
120  ze(i,j,k)=ze(i,j,k)+kmin*d1
121 !
122  xi(i,j,k)=min(xi(i,j,k),1.d0)
123  xi(i,j,k)=max(xi(i,j,k),-1.d0)
124  et(i,j,k)=min(et(i,j,k),1.d0)
125  et(i,j,k)=max(et(i,j,k),-1.d0)
126  ze(i,j,k)=min(ze(i,j,k),1.d0)
127  ze(i,j,k)=max(ze(i,j,k),-1.d0)
128 !
129  call distattach_3d(xi(i,j,k),et(i,j,k),ze(i,j,k),
130  & pneigh,
131  & pnode,a(i,j,k),p,ratio,nterms)
132 c write(*,'(3(1x,i5),4(1x,e11.4))')
133 c & i,j,k,xi(i,j,k),et(i,j,k),ze(i,j,k),a(i,j,k)
134 ! write(*,*) a(i,j,k)
135  endif
136  enddo
137  enddo
138  enddo
139  enddo
140 !
141 ! 2nd run
142 ! initialisation
143 !
144  xiopt=xi(0,0,0)
145  etopt=et(0,0,0)
146  zeopt=ze(0,0,0)
147  do i=-1,1
148  do j=-1,1
149  do k=-1,1
150  xi(i,j,k)=xiopt+i*d2
151  et(i,j,k)=etopt+j*d2
152  ze(i,j,k)=zeopt+k*d2
153  xi(i,j,k)=min(xi(i,j,k),1.d0)
154  xi(i,j,k)=max(xi(i,j,k),-1.d0)
155  et(i,j,k)=min(et(i,j,k),1.d0)
156  et(i,j,k)=max(et(i,j,k),-1.d0)
157  ze(i,j,k)=min(ze(i,j,k),1.d0)
158  ze(i,j,k)=max(ze(i,j,k),-1.d0)
159  call distattach_3d(xi(i,j,k),et(i,j,k),ze(i,j,k),pneigh,
160  & pnode,a(i,j,k),p,ratio,nterms)
161  enddo
162  enddo
163  enddo
164 !
165 ! minimizing the distance from the face to the node
166 !
167  do
168  distmin=a(0,0,0)
169  imin=0
170  jmin=0
171  kmin=0
172  do i=-1,1
173  do j=-1,1
174  do k=-1,1
175  if(a(i,j,k).lt.distmin) then
176  distmin=a(i,j,k)
177  imin=i
178  jmin=j
179  kmin=k
180  endif
181  enddo
182  enddo
183  enddo
184 cc write(*,*) '2 ',xi(0,0,0),et(0,0,0),ze(0,0,0)
185 !
186 ! exit if minimum found
187 !
188  if((imin.eq.0).and.(jmin.eq.0).and.(kmin.eq.0)) exit
189 !
190  do i=-1,1
191  do j=-1,1
192  do k=-1,1
193  aold(i,j,k)=a(i,j,k)
194  xiold(i,j,k)=xi(i,j,k)
195  etold(i,j,k)=et(i,j,k)
196  zeold(i,j,k)=ze(i,j,k)
197  enddo
198  enddo
199  enddo
200 !
201  do i=-1,1
202  do j=-1,1
203  do k=-1,1
204  if((i+imin.ge.-1).and.(i+imin.le.1).and.
205  & (j+jmin.ge.-1).and.(j+jmin.le.1).and.
206  & (k+kmin.ge.-1).and.(k+kmin.le.1)) then
207  a(i,j,k)=aold(i+imin,j+jmin,k+kmin)
208  xi(i,j,k)=xiold(i+imin,j+jmin,k+kmin)
209  et(i,j,k)=etold(i+imin,j+jmin,k+kmin)
210  ze(i,j,k)=zeold(i+imin,j+jmin,k+kmin)
211  else
212  xi(i,j,k)=xi(i,j,k)+imin*d2
213  et(i,j,k)=et(i,j,k)+jmin*d2
214  ze(i,j,k)=ze(i,j,k)+kmin*d2
215 !
216  xi(i,j,k)=min(xi(i,j,k),1.d0)
217  xi(i,j,k)=max(xi(i,j,k),-1.d0)
218  et(i,j,k)=min(et(i,j,k),1.d0)
219  et(i,j,k)=max(et(i,j,k),-1.d0)
220  ze(i,j,k)=min(ze(i,j,k),1.d0)
221  ze(i,j,k)=max(ze(i,j,k),-1.d0)
222 !
223  call distattach_3d(xi(i,j,k),et(i,j,k),ze(i,j,k),
224  & pneigh,
225  & pnode,a(i,j,k),p,ratio,nterms)
226 ! write(*,*) a(i,j,k)
227  endif
228  enddo
229  enddo
230  enddo
231  enddo
232 !
233 ! 3rd run
234 ! initialisation
235 !
236  xiopt=xi(0,0,0)
237  etopt=et(0,0,0)
238  zeopt=ze(0,0,0)
239  do i=-1,1
240  do j=-1,1
241  do k=-1,1
242  xi(i,j,k)=xiopt+i*d3
243  et(i,j,k)=etopt+j*d3
244  ze(i,j,k)=zeopt+k*d3
245  xi(i,j,k)=min(xi(i,j,k),1.d0)
246  xi(i,j,k)=max(xi(i,j,k),-1.d0)
247  et(i,j,k)=min(et(i,j,k),1.d0)
248  et(i,j,k)=max(et(i,j,k),-1.d0)
249  ze(i,j,k)=min(ze(i,j,k),1.d0)
250  ze(i,j,k)=max(ze(i,j,k),-1.d0)
251  call distattach_3d(xi(i,j,k),et(i,j,k),ze(i,j,k),pneigh,
252  & pnode,a(i,j,k),p,ratio,nterms)
253  enddo
254  enddo
255  enddo
256 !
257 ! minimizing the distance from the face to the node
258 !
259  do
260  distmin=a(0,0,0)
261  imin=0
262  jmin=0
263  kmin=0
264  do i=-1,1
265  do j=-1,1
266  do k=-1,1
267  if(a(i,j,k).lt.distmin) then
268  distmin=a(i,j,k)
269  imin=i
270  jmin=j
271  kmin=k
272  endif
273  enddo
274  enddo
275  enddo
276 cc write(*,*) '3 ',xi(0,0,0),et(0,0,0),ze(0,0,0)
277 !
278 ! exit if minimum found
279 !
280  if((imin.eq.0).and.(jmin.eq.0).and.(kmin.eq.0)) exit
281 !
282  do i=-1,1
283  do j=-1,1
284  do k=-1,1
285  aold(i,j,k)=a(i,j,k)
286  xiold(i,j,k)=xi(i,j,k)
287  etold(i,j,k)=et(i,j,k)
288  zeold(i,j,k)=ze(i,j,k)
289  enddo
290  enddo
291  enddo
292 !
293  do i=-1,1
294  do j=-1,1
295  do k=-1,1
296  if((i+imin.ge.-1).and.(i+imin.le.1).and.
297  & (j+jmin.ge.-1).and.(j+jmin.le.1).and.
298  & (k+kmin.ge.-1).and.(k+kmin.le.1)) then
299  a(i,j,k)=aold(i+imin,j+jmin,k+kmin)
300  xi(i,j,k)=xiold(i+imin,j+jmin,k+kmin)
301  et(i,j,k)=etold(i+imin,j+jmin,k+kmin)
302  ze(i,j,k)=zeold(i+imin,j+jmin,k+kmin)
303  else
304  xi(i,j,k)=xi(i,j,k)+imin*d3
305  et(i,j,k)=et(i,j,k)+jmin*d3
306  ze(i,j,k)=ze(i,j,k)+kmin*d3
307 !
308  xi(i,j,k)=min(xi(i,j,k),1.d0)
309  xi(i,j,k)=max(xi(i,j,k),-1.d0)
310  et(i,j,k)=min(et(i,j,k),1.d0)
311  et(i,j,k)=max(et(i,j,k),-1.d0)
312  ze(i,j,k)=min(ze(i,j,k),1.d0)
313  ze(i,j,k)=max(ze(i,j,k),-1.d0)
314 !
315  call distattach_3d(xi(i,j,k),et(i,j,k),ze(i,j,k),
316  & pneigh,
317  & pnode,a(i,j,k),p,ratio,nterms)
318 ! write(*,*) a(i,j,k)
319  endif
320  enddo
321  enddo
322  enddo
323  enddo
324 !
325 ! 4th run
326 ! initialisation
327 !
328  xiopt=xi(0,0,0)
329  etopt=et(0,0,0)
330  zeopt=ze(0,0,0)
331  do i=-1,1
332  do j=-1,1
333  do k=-1,1
334  xi(i,j,k)=xiopt+i*d4
335  et(i,j,k)=etopt+j*d4
336  ze(i,j,k)=zeopt+k*d4
337  xi(i,j,k)=min(xi(i,j,k),1.d0)
338  xi(i,j,k)=max(xi(i,j,k),-1.d0)
339  et(i,j,k)=min(et(i,j,k),1.d0)
340  et(i,j,k)=max(et(i,j,k),-1.d0)
341  ze(i,j,k)=min(ze(i,j,k),1.d0)
342  ze(i,j,k)=max(ze(i,j,k),-1.d0)
343  call distattach_3d(xi(i,j,k),et(i,j,k),ze(i,j,k),pneigh,
344  & pnode,a(i,j,k),p,ratio,nterms)
345  enddo
346  enddo
347  enddo
348 !
349 ! minimizing the distance from the face to the node
350 !
351  do
352  distmin=a(0,0,0)
353  imin=0
354  jmin=0
355  kmin=0
356  do i=-1,1
357  do j=-1,1
358  do k=-1,1
359  if(a(i,j,k).lt.distmin) then
360  distmin=a(i,j,k)
361  imin=i
362  jmin=j
363  kmin=k
364  endif
365  enddo
366  enddo
367  enddo
368 cc write(*,*) '4 ',xi(0,0,0),et(0,0,0),ze(0,0,0)
369 !
370 ! exit if minimum found
371 !
372  if((imin.eq.0).and.(jmin.eq.0).and.(kmin.eq.0)) exit
373 !
374  do i=-1,1
375  do j=-1,1
376  do k=-1,1
377  aold(i,j,k)=a(i,j,k)
378  xiold(i,j,k)=xi(i,j,k)
379  etold(i,j,k)=et(i,j,k)
380  zeold(i,j,k)=ze(i,j,k)
381  enddo
382  enddo
383  enddo
384 !
385  do i=-1,1
386  do j=-1,1
387  do k=-1,1
388  if((i+imin.ge.-1).and.(i+imin.le.1).and.
389  & (j+jmin.ge.-1).and.(j+jmin.le.1).and.
390  & (k+kmin.ge.-1).and.(k+kmin.le.1)) then
391  a(i,j,k)=aold(i+imin,j+jmin,k+kmin)
392  xi(i,j,k)=xiold(i+imin,j+jmin,k+kmin)
393  et(i,j,k)=etold(i+imin,j+jmin,k+kmin)
394  ze(i,j,k)=zeold(i+imin,j+jmin,k+kmin)
395  else
396  xi(i,j,k)=xi(i,j,k)+imin*d4
397  et(i,j,k)=et(i,j,k)+jmin*d4
398  ze(i,j,k)=ze(i,j,k)+kmin*d4
399 !
400  xi(i,j,k)=min(xi(i,j,k),1.d0)
401  xi(i,j,k)=max(xi(i,j,k),-1.d0)
402  et(i,j,k)=min(et(i,j,k),1.d0)
403  et(i,j,k)=max(et(i,j,k),-1.d0)
404  ze(i,j,k)=min(ze(i,j,k),1.d0)
405  ze(i,j,k)=max(ze(i,j,k),-1.d0)
406 !
407  call distattach_3d(xi(i,j,k),et(i,j,k),ze(i,j,k),
408  & pneigh,
409  & pnode,a(i,j,k),p,ratio,nterms)
410 ! write(*,*) a(i,j,k)
411  endif
412  enddo
413  enddo
414  enddo
415  enddo
416 !
417  call distattach_3d(xi(0,0,0),et(0,0,0),ze(0,0,0),pneigh,pnode,
418  & a(0,0,0),p,ratio,nterms)
419 !
420  do i=1,3
421  pnode(i)=p(i)
422  enddo
423 !
424  dist=a(0,0,0)
425 !
426  if((nterms.eq.20).or.(nterms.eq.8)) then
427  xil=xi(0,0,0)
428  etl=et(0,0,0)
429  zel=ze(0,0,0)
430  elseif((nterms.eq.4).or.(nterms.eq.10)) then
431  xil=(xi(0,0,0)+1.d0)/2.d0
432  etl=(et(0,0,0)+1.d0)/2.d0
433  zel=(ze(0,0,0)+1.d0)/2.d0
434  dx(1)=xil
435  dx(2)=etl
436  dx(3)=zel
437  call dsort(dx,iy,n,kflag)
438  if(dx(3).gt.1.d-30) then
439  al=dx(3)/(xil+etl+zel)
440  xil=al*xil
441  etl=al*etl
442  zel=al*zel
443  endif
444 c if(xil+etl+zel.gt.1.d0) then
445 c dummy=2.d0*(1.d0-xil-etl-zel)/3.d0
446 c xil=dummy+xil
447 c etl=dummy+etl
448 c zel=dummy+zel
449 c endif
450  elseif((nterms.eq.6).or.(nterms.eq.15)) then
451  xil=(xi(0,0,0)+1.d0)/2.d0
452  etl=(et(0,0,0)+1.d0)/2.d0
453  if(xil+etl.gt.1.d0) then
454  dummy=xil
455  xil=1.d0-etl
456  etl=1.d0-dummy
457  endif
458  zel=ze(0,0,0)
459  endif
460 !
461  return
#define max(a, b)
Definition: cascade.c:32
#define min(a, b)
Definition: cascade.c:31
static double * dist
Definition: radflowload.c:42
subroutine dsort(dx, iy, n, kflag)
Definition: dsort.f:6
subroutine distattach_3d(xig, etg, zeg, pneigh, pnode, a, p, ratio, nterms)
Definition: distattach_3d.f:21
Hosted by OpenAircraft.com, (Michigan UAV, LLC)