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

Go to the source code of this file.

Functions/Subroutines

subroutine map3dtolayer (yn, ipkon, kon, lakon, nfield, ne, co, ielmat, mi)
 

Function/Subroutine Documentation

◆ map3dtolayer()

subroutine map3dtolayer ( real*8, dimension(nfield,*)  yn,
integer, dimension(*)  ipkon,
integer, dimension(*)  kon,
character*8, dimension(*)  lakon,
integer  nfield,
integer  ne,
real*8, dimension(3,*)  co,
integer, dimension(mi(3),*)  ielmat,
integer, dimension(*)  mi 
)
21 !
22 ! interpolates 3d field nodal values to nodal values in the
23 ! layers of composite materials
24 !
25  implicit none
26 !
27  character*8 lakon(*)
28 !
29  integer ipkon(*),kon(*),ne,nfield,i,j,k,l,mi(*),nelem,
30  & ielmat(mi(3),*),nlayer,iflag,nbot20(8),nmid20(8),ntop20(8),
31  & nbot15(6),nmid15(6),ntop15(6),ibot,itop,nope,nopes,nopeexp,
32  & indexe,konl(20),imid,node,nodebot,nodetop
33 !
34  real*8 yn(nfield,*),shp(4,20),xsj(3),co(3,*),xl(3,20),
35  & xi,et,ze,dd,dt,xi20(8),et20(8),xi15(6),et15(6)
36 !
37  data iflag /1/
38 !
39  data nbot20 /1,2,3,4,9,10,11,12/
40  data nmid20 /17,18,19,20,0,0,0,0/
41  data ntop20 /5,6,7,8,13,14,15,16/
42 !
43  data nbot15 /1,2,3,7,8,9/
44  data nmid15 /13,14,15,0,0,0/
45  data ntop15 /4,5,6,10,11,12/
46 !
47  data xi20 /-1.d0,1.d0,1.d0,-1.d0,0.d0,1.d0,0.d0,-1.d0/
48  data et20 /-1.d0,-1.d0,1.d0,1.d0,-1.d0,0.d0,1.d0,0.d0/
49 !
50  data xi15 /0.d0,1.d0,0.d0,0.5d0,0.5d0,0.d0/
51  data et15 /0.d0,0.d0,1.d0,0.d0,0.5d0,0.5d0/
52 !
53  include "gauss.f"
54 !
55  do nelem=1,ne
56  if(lakon(nelem)(7:8).eq.'LC') then
57 !
58 ! composite materials
59 !
60 ! determining the number of layers
61 !
62  nlayer=0
63  do k=1,mi(3)
64  if(ielmat(k,nelem).ne.0) then
65  nlayer=nlayer+1
66  endif
67  enddo
68 !
69  indexe=ipkon(nelem)
70 !
71  if(lakon(nelem)(4:5).eq.'20') then
72  nope=20
73  nopes=8
74  nopeexp=28
75  elseif(lakon(nelem)(4:5).eq.'15') then
76  nope=15
77  nopes=6
78  nopeexp=21
79  endif
80 !
81  do i=1,nope
82  konl(i)=kon(indexe+i)
83  do j=1,3
84  xl(j,i)=co(j,konl(i))
85  enddo
86  enddo
87 !
88  do i=1,nopes
89  if(lakon(nelem)(4:5).eq.'20') then
90  ibot=nbot20(i)
91  imid=nmid20(i)
92  itop=ntop20(i)
93  xi=xi20(i)
94  et=et20(i)
95  else
96  ibot=nbot15(i)
97  imid=nmid15(i)
98  itop=ntop15(i)
99  xi=xi15(i)
100  et=et15(i)
101  endif
102 !
103  nodebot=konl(ibot)
104  nodetop=konl(itop)
105  dd=sqrt((co(1,nodebot)-co(1,nodetop))**2+
106  & (co(2,nodebot)-co(2,nodetop))**2+
107  & (co(3,nodebot)-co(3,nodetop))**2)
108  do j=0,nlayer-1
109 !
110 ! bottom node
111 !
112  node=kon(indexe+nopeexp+j*nope+ibot)
113  dt=sqrt((co(1,nodebot)-co(1,node))**2+
114  & (co(2,nodebot)-co(2,node))**2+
115  & (co(3,nodebot)-co(3,node))**2)
116  ze=2.d0*dt/dd-1.d0
117 c write(*,*) 'map3dtolayer',node,xi,et,ze
118 !
119 ! determining the value of the shape functions
120 !
121  if(lakon(nelem)(4:5).eq.'20') then
122  call shape20h(xi,et,ze,xl,xsj,shp,iflag)
123  elseif(lakon(nelem)(4:5).eq.'15') then
124  call shape15w(xi,et,ze,xl,xsj,shp,iflag)
125  endif
126 !
127  do k=1,nfield
128  yn(k,node)=0.d0
129  do l=1,nope
130  yn(k,node)=yn(k,node)+
131  & shp(4,l)*yn(k,konl(l))
132 c write(*,*) 'xi',xi,et,ze,node
133 c write(*,*) l,shp(4,l),yn(k,konl(l))
134  enddo
135  enddo
136 !
137 ! top node
138 !
139  node=kon(indexe+nopeexp+j*nope+itop)
140  dt=sqrt((co(1,nodebot)-co(1,node))**2+
141  & (co(2,nodebot)-co(2,node))**2+
142  & (co(3,nodebot)-co(3,node))**2)
143  ze=2.d0*dt/dd-1.d0
144 c write(*,*) 'map3dtolayer',node,xi,et,ze
145 !
146 ! determining the value of the shape functions
147 !
148  if(lakon(nelem)(4:5).eq.'20') then
149  call shape20h(xi,et,ze,xl,xsj,shp,iflag)
150  elseif(lakon(nelem)(4:5).eq.'15') then
151  call shape15w(xi,et,ze,xl,xsj,shp,iflag)
152  endif
153 !
154  do k=1,nfield
155  yn(k,node)=0.d0
156  do l=1,nope
157  yn(k,node)=yn(k,node)+
158  & shp(4,l)*yn(k,konl(l))
159  enddo
160  enddo
161 !
162 ! middle node, if any
163 !
164  if(i.gt.nopes/2) cycle
165 !
166  node=kon(indexe+nopeexp+j*nope+imid)
167  dt=sqrt((co(1,nodebot)-co(1,node))**2+
168  & (co(2,nodebot)-co(2,node))**2+
169  & (co(3,nodebot)-co(3,node))**2)
170  ze=2.d0*dt/dd-1.d0
171 c write(*,*) 'map3dtolayer',node,xi,et,ze
172 !
173 ! determining the value of the shape functions
174 !
175  if(lakon(nelem)(4:5).eq.'20') then
176  call shape20h(xi,et,ze,xl,xsj,shp,iflag)
177  elseif(lakon(nelem)(4:5).eq.'15') then
178  call shape15w(xi,et,ze,xl,xsj,shp,iflag)
179  endif
180 !
181  do k=1,nfield
182  yn(k,node)=0.d0
183  do l=1,nope
184  yn(k,node)=yn(k,node)+
185  & shp(4,l)*yn(k,konl(l))
186  enddo
187  enddo
188  enddo
189  enddo
190  endif
191  enddo
192 !
193  return
subroutine shape15w(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape15w.f:20
subroutine shape20h(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape20h.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)