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

Go to the source code of this file.

Functions/Subroutines

subroutine gen3dfrom2d (i, kon, ipkon, lakon, ne, iponor, xnor, knor, thicke, offset, ntrans, inotr, trab, ikboun, ilboun, nboun, nboun_, nodeboun, ndirboun, xboun, iamboun, typeboun, ipompc, nodempc, coefmpc, nmpc, nmpc_, mpcfree, ikmpc, ilmpc, labmpc, nk, nk_, co, rig, nmethod, iperturb, ithermal, mi, nam, icomposite, ielmat, vold, orname, orab, norien, norien_, ielorien)
 

Function/Subroutine Documentation

◆ gen3dfrom2d()

subroutine gen3dfrom2d ( integer  i,
integer, dimension(*)  kon,
integer, dimension(*)  ipkon,
character*8, dimension(*)  lakon,
integer  ne,
integer, dimension(2,*)  iponor,
real*8, dimension(*)  xnor,
integer, dimension(*)  knor,
real*8, dimension(mi(3),*)  thicke,
real*8, dimension(2,*)  offset,
integer  ntrans,
integer, dimension(2,*)  inotr,
real*8, dimension(7,*)  trab,
integer, dimension(*)  ikboun,
integer, dimension(*)  ilboun,
integer  nboun,
integer  nboun_,
integer, dimension(*)  nodeboun,
integer, dimension(*)  ndirboun,
real*8, dimension(*)  xboun,
integer, dimension(*)  iamboun,
character*1, dimension(*)  typeboun,
integer, dimension(*)  ipompc,
integer, dimension(3,*)  nodempc,
real*8, dimension(*)  coefmpc,
integer  nmpc,
integer  nmpc_,
integer  mpcfree,
integer, dimension(*)  ikmpc,
integer, dimension(*)  ilmpc,
character*20, dimension(*)  labmpc,
integer  nk,
integer  nk_,
real*8, dimension(3,*)  co,
integer, dimension(*)  rig,
integer  nmethod,
integer  iperturb,
integer, dimension(2)  ithermal,
integer, dimension(*)  mi,
integer  nam,
integer  icomposite,
integer, dimension(mi(3),*)  ielmat,
real*8, dimension(0:mi(2),*)  vold,
character*80, dimension(*)  orname,
real*8, dimension(7,*)  orab,
integer  norien,
integer  norien_,
integer, dimension(mi(3),*)  ielorien 
)
25 !
26 ! expands 2d element i into a 3d element
27 !
28 ! generates additional MPC's for plane stress, plane strain and
29 ! axisymmetric elements
30 !
31  implicit none
32 !
33  logical axial,fixed,quadratic
34 !
35  character*1 type,typeboun(*)
36  character*8 lakon(*)
37  character*20 labmpc(*),label
38  character*80 orname(*)
39 !
40  integer kon(*),ipkon(*),ne,iponor(2,*),knor(*),ntrans,
41  & inotr(2,*),nodefixz,norien,norien_,norieninput,iflag,
42  & ikboun(*),ilboun(*),nboun,nboun_,nodeboun(*),ndirboun(*),
43  & iamboun(*),nam,ipompc(*),nodempc(3,*),nmpc,nmpc_,mpcfree,
44  & ikmpc(*),ilmpc(*),nk,nk_,i,rig(*),nmethod,iperturb,ishift,
45  & indexe,j,nodel(8),indexx,indexk,k,nedge,nodes(3,8),nodec(3,8),
46  & iamplitude,l,newnode,idir,idof,id,m,mpcfreenew,node,ithermal(2),
47  & jmin,jmax,idummy,mi(*),indexc,indexl,icomposite,ielmat(mi(3),*),
48  & nope,neworien(0:norien),iorien,ipos,ielorien(mi(3),*)
49 !
50  real*8 xnor(*),thicke(mi(3),*),offset(2,*),trab(7,*),xboun(*),
51  & coefmpc(*),co(3,*),thicks(8),xnors(3,8),dc,ds,val,orab(7,*),
52  & x,y,thickness(8),vold(0:mi(2),*),xi,et,xl(3,8),xno(3),xs(3,7),
53  & shp(7,8),p(3),a(3,3),dd,e1(3),e2(3)
54 !
55  norieninput=norien
56  iflag=2
57 !
58  label=' '
59  fixed=.false.
60 !
61  if(lakon(i)(1:2).eq.'CA') then
62  axial=.true.
63  else
64  axial=.false.
65  endif
66 !
67  indexe=ipkon(i)
68 !
69 ! localizing the nodes, thicknesses and normals for the
70 ! 2-D element
71 !
72  if((lakon(i)(2:2).eq.'3').or.
73  & (lakon(i)(4:4).eq.'3')) then
74  quadratic=.false.
75  nedge=3
76  nope=3
77  ishift=6
78  elseif((lakon(i)(2:2).eq.'4').or.
79  & (lakon(i)(4:4).eq.'4')) then
80  quadratic=.false.
81  nedge=4
82  nope=4
83 !
84 ! CPE4R,CPS4R, CAX4R and S4R elements are expanded into C3D8R
85 !
86  if((lakon(i)(3:3).eq.'R').or.
87  & (lakon(i)(5:5).eq.'R')) then
88  ishift=8
89 !
90 ! S4 elements are expanded into C3D8I
91 !
92  elseif((lakon(i)(1:1).eq.'S').or.(lakon(i)(1:1).eq.'M')) then
93  ishift=11
94 !
95 ! CPE4, CPS4 and CAX4 elements are expanded into C3D8
96 !
97  else
98  ishift=8
99  endif
100  elseif((lakon(i)(2:2).eq.'6').or.
101  & (lakon(i)(4:4).eq.'6')) then
102  quadratic=.true.
103  nedge=3
104  nope=6
105  ishift=15
106  elseif((lakon(i)(2:2).eq.'8').or.
107  & (lakon(i)(4:4).eq.'8')) then
108  quadratic=.true.
109  nedge=4
110  nope=8
111  ishift=20
112  endif
113 !
114 c do j=1,2*nedge
115  do j=1,nope
116  nodel(j)=kon(indexe+j)
117  kon(indexe+ishift+j)=nodel(j)
118  indexk=iponor(2,indexe+j)
119  thicks(j)=0.d0
120  do k=1,mi(3)
121  thicks(j)=thicks(j)+thicke(k,indexe+j)
122  enddo
123  do k=1,3
124  nodes(k,j)=knor(indexk+k)
125  enddo
126  enddo
127 !
128 ! generating the 3-D element topology for shell and plane
129 ! stress/strain elements
130 !
131  if(lakon(i)(1:2).ne.'CA') then
132 c do j=1,2*nedge
133  do j=1,nope
134  indexx=iponor(1,indexe+j)
135  do k=1,3
136  xnors(k,j)=xnor(indexx+k)
137  enddo
138  if(ntrans.gt.0) then
139  do k=1,3
140  inotr(1,nodes(k,j))=inotr(1,nodel(j))
141  enddo
142  endif
143  enddo
144 !
145 ! vertex nodes
146 !
147  do k=1,nedge
148  kon(indexe+k)=nodes(1,k)
149 !
150  do j=1,3
151  co(j,nodes(1,k))=co(j,nodel(k))
152  & -thicks(k)*xnors(j,k)*(.5d0+offset(1,i))
153  enddo
154  enddo
155  do k=1,nedge
156  kon(indexe+nedge+k)=nodes(3,k)
157  do j=1,3
158  co(j,nodes(3,k))=co(j,nodel(k))
159  & +thicks(k)*xnors(j,k)*(.5d0-offset(1,i))
160  enddo
161  enddo
162 !
163 ! middle nodes for quadratic expansion
164 !
165  if(quadratic) then
166  do k=nedge+1,2*nedge
167  kon(indexe+nedge+k)=nodes(1,k)
168  do j=1,3
169  co(j,nodes(1,k))=co(j,nodel(k))
170  & -thicks(k)*xnors(j,k)*(.5d0+offset(1,i))
171  enddo
172  enddo
173  do k=nedge+1,2*nedge
174  kon(indexe+2*nedge+k)=nodes(3,k)
175  do j=1,3
176  co(j,nodes(3,k))=co(j,nodel(k))
177  & +thicks(k)*xnors(j,k)*(.5d0-offset(1,i))
178  enddo
179  enddo
180  do k=1,nedge
181  kon(indexe+4*nedge+k)=nodes(2,k)
182  do j=1,3
183  co(j,nodes(2,k))=co(j,nodel(k))
184  & -thicks(k)*xnors(j,k)*offset(1,i)
185  enddo
186  enddo
187 !
188 ! generating coordinates for the expanded nodes which
189 ! are not used by the C3D20(R) element (needed for the
190 ! determination of the knot dimension)
191 !
192  do k=nedge+1,2*nedge
193  do j=1,3
194  co(j,nodes(2,k))=co(j,nodel(k))
195  & -thicks(k)*xnors(j,k)*offset(1,i)
196  enddo
197  enddo
198  else
199 !
200 ! for linear elements the coordinates of the nodes
201 ! in the middle may be needed in case of knots (to
202 ! determine the coefficients of the equations)
203 !
204  do k=1,nedge
205  do j=1,3
206  co(j,nodes(2,k))=co(j,nodel(k))
207  & -thicks(k)*xnors(j,k)*offset(1,i)
208  enddo
209  enddo
210  endif
211 !
212 ! generating the layer geometry for composite shells
213 ! only for S8R elements
214 !
215  if(lakon(i)(8:8).eq.'C') then
216  do k=1,2*nedge
217  thickness(k)=0.d0
218  enddo
219  indexc=indexe+2*nedge
220  indexl=0
221  do l=1,mi(3)
222  if(ielmat(l,i).eq.0) exit
223  indexc=indexc+ishift
224  indexl=indexl+3
225 !
226  do j=1,2*nedge
227  do k=1,3
228  nodec(k,j)=knor(iponor(2,indexe+j)+indexl+k)
229  enddo
230  enddo
231 !
232 ! nodes 1 to 4 (vertex nodes face 1)
233 !
234  do k=1,nedge
235  kon(indexc+k)=nodec(1,k)
236  do j=1,3
237  co(j,nodec(1,k))=co(j,nodel(k))+
238  & (thickness(k)-thicks(k)*(.5d0+offset(1,i)))
239  & *xnors(j,k)
240  enddo
241  enddo
242 !
243 ! nodes 5 to 8 (vertex nodes face 2)
244 !
245  do k=1,nedge
246  kon(indexc+nedge+k)=nodec(3,k)
247  do j=1,3
248  co(j,nodec(3,k))=co(j,nodel(k))+
249  & (thickness(k)+thicke(l,indexe+k)
250  & -thicks(k)*(.5d0+offset(1,i)))
251  & *xnors(j,k)
252  enddo
253  enddo
254 !
255 ! nodes 9 to 12 (middle nodes face 1)
256 !
257  do k=nedge+1,2*nedge
258  kon(indexc+nedge+k)=nodec(1,k)
259  do j=1,3
260  co(j,nodec(1,k))=co(j,nodel(k))+
261  & (thickness(k)-thicks(k)*(.5d0+offset(1,i)))
262  & *xnors(j,k)
263  enddo
264  enddo
265 !
266 ! nodes 13 to 16 (middle nodes face 2)
267 !
268  do k=nedge+1,2*nedge
269  kon(indexc+2*nedge+k)=
270  & nodec(3,k)
271  do j=1,3
272  co(j,nodec(3,k))=co(j,nodel(k))+
273  & (thickness(k)+thicke(l,indexe+k)
274  & -thicks(k)*(.5d0+offset(1,i)))
275  & *xnors(j,k)
276  enddo
277  enddo
278 !
279 ! nodes 17 to 20 (middle nodes in between face 1 and 2)
280 !
281  do k=1,nedge
282  kon(indexc+4*nedge+k)=
283  & nodec(2,k)
284  do j=1,3
285  co(j,nodec(2,k))=co(j,nodel(k))+
286  & (thickness(k)+thicke(l,indexe+k)/2.d0
287  & -thicks(k)*(.5d0+offset(1,i)))
288  & *xnors(j,k)
289  enddo
290  enddo
291 !
292  do k=1,2*nedge
293  thickness(k)=thickness(k)+thicke(l,indexe+k)
294  enddo
295  enddo
296  endif
297  else
298 !
299 ! generating the 3-D element topology for axisymmetric elements
300 !
301  dc=dcos(thicks(1)/2.d0)
302  ds=dsin(thicks(1)/2.d0)
303  do j=1,nedge
304  indexk=iponor(2,indexe+j)
305  x=co(1,nodel(j))
306  y=co(2,nodel(j))
307 !
308  node=knor(indexk+1)
309  co(1,node)=x*dc
310  co(2,node)=y
311  co(3,node)=-x*ds
312  kon(indexe+j)=node
313 !
314  if(quadratic) then
315  node=knor(indexk+2)
316  co(1,node)=x
317  co(2,node)=y
318  co(3,node)=0.d0
319  kon(indexe+4*nedge+j)=node
320  else
321 !
322 ! for linear elements the coordinates of the nodes
323 ! in the middle may be needed in case of knots (to
324 ! determine the coefficients of the equations)
325 !
326  node=knor(indexk+2)
327  co(1,node)=x
328  co(2,node)=y
329  co(3,node)=0.d0
330  endif
331 !
332  node=knor(indexk+3)
333  co(1,node)=x*dc
334  co(2,node)=y
335  co(3,node)=x*ds
336  kon(indexe+nedge+j)=node
337  enddo
338 !
339  if(quadratic) then
340  do j=nedge+1,2*nedge
341  indexk=iponor(2,indexe+j)
342  x=co(1,nodel(j))
343  y=co(2,nodel(j))
344 !
345  node=knor(indexk+1)
346  co(1,node)=x*dc
347  co(2,node)=y
348  co(3,node)=-x*ds
349  kon(indexe+nedge+j)=node
350 !
351  node=knor(indexk+3)
352  co(1,node)=x*dc
353  co(2,node)=y
354  co(3,node)=x*ds
355  kon(indexe+2*nedge+j)=node
356  enddo
357  endif
358  endif
359 !
360 ! additional SPC's due to plane strain/plane stress/axisymmetric
361 ! conditions
362 !
363  do j=1,nope
364  if((lakon(i)(1:1).ne.'S').and.(lakon(i)(1:1).ne.'M')) then
365 !
366 ! fixing the middle plane
367 !
368 c if(rig(nodel(j)).gt.0) cycle
369  if(rig(nodel(j)).gt.0) then
370  nodefixz=nodel(j)
371  else
372  nodefixz=nodes(2,j)
373  endif
374 !
375  if(ithermal(2).ne.2) then
376  val=0.d0
377  k=3
378  if(nam.gt.0) iamplitude=0
379  type='M'
380 c call bounadd(nodes(2,j),k,k,val,nodeboun,
381  call bounadd(nodefixz,k,k,val,nodeboun,
382  & ndirboun,xboun,nboun,nboun_,iamboun,
383  & iamplitude,nam,ipompc,nodempc,coefmpc,
384  & nmpc,nmpc_,mpcfree,inotr,trab,ntrans,
385  & ikboun,ilboun,ikmpc,ilmpc,co,nk,nk_,
386  & labmpc,type,typeboun,nmethod,iperturb,
387  & fixed,vold,idummy,mi,label)
388  endif
389 !
390 ! specifying that the side planes do the same
391 ! as the middle plane (in all directions for
392 ! plane strain and axisymmetric elements, in the
393 ! plane for plane stress elements)
394 !
395  if(ithermal(2).le.1) then
396  jmin=1
397  jmax=3
398  elseif(ithermal(2).eq.2) then
399  jmin=0
400  jmax=0
401  else
402  jmin=0
403  jmax=3
404  endif
405 !
406  do l=1,3,2
407  newnode=nodes(l,j)
408  do idir=jmin,jmax
409  if((idir.eq.3).and.(lakon(i)(1:3).eq.'CPS')) then
410 !
411 ! for linear elements and plane stress it is not
412 ! sufficient to set the z-displacement of the
413 ! midplane nodes to zero, since these are not used.
414 ! Instead, the sum of the z-displacement of the
415 ! boundary plane nodes are set to zero:
416 ! w_1+w_3=0
417 !
418  if(((lakon(i)(4:4).eq.'3').or.
419  & (lakon(i)(4:4).eq.'4')).and.
420  & (l.eq.1)) then
421  idof=8*(newnode-1)+3
422  call nident(ikmpc,idof,nmpc,id)
423  if((id.le.0).or.(ikmpc(id).ne.idof)) then
424  nmpc=nmpc+1
425  if(nmpc.gt.nmpc_) then
426  write(*,*)
427  & '*ERROR in gen3dfrom2d: increase nmpc_'
428  call exit(201)
429  endif
430  labmpc(nmpc)=' '
431  ipompc(nmpc)=mpcfree
432  do m=nmpc,id+2,-1
433  ikmpc(m)=ikmpc(m-1)
434  ilmpc(m)=ilmpc(m-1)
435  enddo
436  ikmpc(id+1)=idof
437  ilmpc(id+1)=nmpc
438  nodempc(1,mpcfree)=newnode
439  nodempc(2,mpcfree)=3
440  coefmpc(mpcfree)=1.d0
441  mpcfree=nodempc(3,mpcfree)
442  if(mpcfree.eq.0) then
443  write(*,*)
444  & '*ERROR in gen3dfrom2d: increase memmpc_'
445  call exit(201)
446  endif
447  nodempc(1,mpcfree)=nodes(3,j)
448  nodempc(2,mpcfree)=3
449  coefmpc(mpcfree)=1.d0
450  mpcfreenew=nodempc(3,mpcfree)
451  if(mpcfreenew.eq.0) then
452  write(*,*)
453  & '*ERROR in gen3dfrom2d: increase memmpc_'
454  call exit(201)
455  endif
456  nodempc(3,mpcfree)=0
457  mpcfree=mpcfreenew
458  endif
459  else
460  cycle
461  endif
462  endif
463  idof=8*(newnode-1)+idir
464  call nident(ikmpc,idof,nmpc,id)
465  if((id.le.0).or.(ikmpc(id).ne.idof)) then
466  nmpc=nmpc+1
467  if(nmpc.gt.nmpc_) then
468  write(*,*)
469  & '*ERROR in gen3dfrom2d: increase nmpc_'
470  call exit(201)
471  endif
472  labmpc(nmpc)=' '
473  ipompc(nmpc)=mpcfree
474  do m=nmpc,id+2,-1
475  ikmpc(m)=ikmpc(m-1)
476  ilmpc(m)=ilmpc(m-1)
477  enddo
478  ikmpc(id+1)=idof
479  ilmpc(id+1)=nmpc
480  nodempc(1,mpcfree)=newnode
481  nodempc(2,mpcfree)=idir
482  coefmpc(mpcfree)=1.d0
483  mpcfree=nodempc(3,mpcfree)
484  if(mpcfree.eq.0) then
485  write(*,*)
486  & '*ERROR in gen3dfrom2d: increase memmpc_'
487  call exit(201)
488  endif
489  nodempc(1,mpcfree)=nodes(2,j)
490  if((lakon(i)(2:2).eq.'A').and.(idir.eq.3))
491  & then
492  nodempc(2,mpcfree)=1
493  else
494  nodempc(2,mpcfree)=idir
495  endif
496  if(lakon(i)(2:2).eq.'A') then
497  if(idir.eq.1) then
498  coefmpc(mpcfree)=-dc
499  elseif(idir.eq.3) then
500  if(l.eq.1) then
501  coefmpc(mpcfree)=ds
502  else
503  coefmpc(mpcfree)=-ds
504  endif
505  else
506  coefmpc(mpcfree)=-1.d0
507  endif
508  else
509  coefmpc(mpcfree)=-1.d0
510  endif
511  mpcfreenew=nodempc(3,mpcfree)
512  if(mpcfreenew.eq.0) then
513  write(*,*)
514  & '*ERROR in gen3dfrom2d: increase memmpc_'
515  call exit(201)
516  endif
517  nodempc(3,mpcfree)=0
518  mpcfree=mpcfreenew
519  endif
520  enddo
521  enddo
522  endif
523  enddo
524 !
525 ! create a local coordinate system for each layer in a shell
526 ! element
527 !
528  if(lakon(i)(1:1).eq.'S') then
529 !
530 ! initialize neworien
531 !
532  do j=0,norieninput
533  neworien(j)=0
534  enddo
535 !
536 ! loop over the layers
537 !
538  do j=1,mi(3)
539 !
540 ! check whether layer exists (i.e. whether a material
541 ! was assigned to the layer)
542 !
543  if(ielmat(j,i).le.0) exit
544 !
545  iorien=ielorien(j,i)
546 !
547 ! check whether this orientation already occurred within
548 ! this element (i.e. for a previous layer)
549 !
550  if(neworien(iorien).gt.0) then
551  ielorien(j,i)=neworien(iorien)
552  cycle
553  endif
554 !
555 ! calculate the coordinates of the point at
556 ! (xi,et)=(0,0) and the normal at that point
557 !
558  xi=0.d0
559  et=0.d0
560 !
561  if(lakon(i)(2:2).eq.'3') then
562  nope=3
563  elseif(lakon(i)(2:2).eq.'4') then
564  nope=4
565  elseif(lakon(i)(2:2).eq.'6') then
566  nope=6
567  elseif(lakon(i)(2:2).eq.'8') then
568  nope=8
569  else
570  exit
571  endif
572 !
573  do k=1,nope
574  node=kon(indexe+ishift+k)
575  do l=1,3
576  xl(l,k)=co(l,node)
577  enddo
578  enddo
579 !
580  if(nope.eq.3) then
581  call shape3tri(xi,et,xl,xno,xs,shp,iflag)
582  elseif(nope.eq.4) then
583  call shape4q(xi,et,xl,xno,xs,shp,iflag)
584  elseif(nope.eq.6) then
585  call shape6tri(xi,et,xl,xno,xs,shp,iflag)
586  else
587  call shape8q(xi,et,xl,xno,xs,shp,iflag)
588  endif
589 !
590  dd=dsqrt(xno(1)*xno(1)+xno(2)*xno(2)+xno(3)*xno(3))
591  do l=1,3
592  xno(l)=xno(l)/dd
593  enddo
594 c write(*,*) 'gen3dfrom2d xnor',(xnor(l),l=1,3)
595 !
596 ! coordinates of the point at (xi,et)=(0.,0.)
597 !
598  do l=1,3
599  p(l)=0.d0
600  do k=1,nope
601  p(l)=p(l)+shp(4,k)*xl(l,k)
602  enddo
603  enddo
604 !
605  if(iorien.eq.0) then
606 !
607 ! unit matrix
608 !
609  do k=1,3
610  do l=1,3
611  a(k,l)=0.d0
612  enddo
613  a(k,k)=1.d0
614  enddo
615  else
616  call transformatrix(orab(1,iorien),p,a)
617  endif
618 !
619 ! dd=e1.n
620 !
621  dd=a(1,1)*xno(1)+a(2,1)*xno(2)+a(3,1)*xno(3)
622 !
623 ! check whether e1 is within 0.1° from n
624 ! e1 is the local 1-direction in the shell plane
625 !
626  if(dabs(dd).gt.0.999999999536d0) then
627 !
628 ! project the z-axis
629 !
630  dd=a(1,3)*xno(1)+a(2,3)*xno(2)+a(3,3)*xno(3)
631  do l=1,3
632  e1(l)=a(l,3)-dd*xno(l)
633  enddo
634  else
635 !
636 ! project the x-axis
637 !
638  do l=1,3
639  e1(l)=a(l,1)-dd*xno(l)
640  enddo
641  endif
642 !
643  dd=dsqrt(e1(1)*e1(1)+e1(2)*e1(2)+e1(3)*e1(3))
644  do l=1,3
645  e1(l)=e1(l)/dd
646  enddo
647 !
648 ! e2=nxe1
649 !
650  e2(1)=xno(2)*e1(3)-e1(2)*xno(3)
651  e2(2)=xno(3)*e1(1)-e1(3)*xno(1)
652  e2(3)=xno(1)*e1(2)-e1(1)*xno(2)
653 !
654  norien=norien+1
655  if(norien.gt.norien_) then
656  write(*,*) '*ERROR in gen3dfrom2d: increase norien_'
657  call exit(201)
658  endif
659 !
660 ! creating the name of the new orientation =
661 ! old name + _shell_ + element number (10 digits;
662 ! this corresponds to the maximum number of digits
663 ! for element numbers in the frd-format)
664 !
665  if(iorien.gt.0) then
666  ipos=index(orname(iorien),' ')
667  orname(norien)(1:ipos-1)=orname(iorien)(1:ipos-1)
668  else
669  ipos=1
670  endif
671  orname(norien)(ipos:ipos+6)='_shell_'
672  write(orname(norien)(ipos+7:ipos+16),'(i10.10)') i
673  do l=ipos+17,80
674  orname(norien)(l:l)=' '
675  enddo
676 !
677  do l=1,3
678  orab(l,norien)=e1(l)
679  orab(l+3,norien)=e2(l)
680  enddo
681  orab(7,norien)=1.d0
682 c write(*,*) 'gen3dfrom2d o',(orab(l,iorien),l=1,7)
683 c write(*,*) 'gen3dfrom2d n',(orab(l,norien),l=1,7)
684 !
685  neworien(iorien)=norien
686  ielorien(j,i)=norien
687 !
688  enddo
689 !
690  endif
691 !
692  return
subroutine shape8q(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape8q.f:20
subroutine shape3tri(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape3tri.f:20
subroutine transformatrix(xab, p, a)
Definition: transformatrix.f:20
subroutine shape4q(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape4q.f:20
subroutine nident(x, px, n, id)
Definition: nident.f:26
subroutine nodes(inpc, textpart, co, nk, nk_, set, istartset, iendset, ialset, nset, nset_, nalset, nalset_, istep, istat, n, iline, ipol, inl, ipoinp, inp, ipoinpc)
Definition: nodes.f:22
subroutine bounadd(node, is, ie, val, nodeboun, ndirboun, xboun, nboun, nboun_, iamboun, iamplitude, nam, ipompc, nodempc, coefmpc, nmpc, nmpc_, mpcfree, inotr, trab, ntrans, ikboun, ilboun, ikmpc, ilmpc, co, nk, nk_, labmpc, type, typeboun, nmethod, iperturb, fixed, vold, nodetrue, mi, label)
Definition: bounadd.f:24
subroutine shape6tri(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape6tri.f:20
subroutine thickness(dgdx, nobject, nodedesiboun, ndesiboun, objectset, xo, yo, zo, x, y, z, nx, ny, nz, co, ifree, ndesia, ndesib, iobject, ndesi, dgdxglob, nk)
Definition: thickness.f:22
Hosted by OpenAircraft.com, (Michigan UAV, LLC)