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

Go to the source code of this file.

Functions/Subroutines

subroutine map3dto1d2d_v (yn, ipkon, inum, kon, lakon, nfield, nk, ne, nactdof)
 

Function/Subroutine Documentation

◆ map3dto1d2d_v()

subroutine map3dto1d2d_v ( real*8, dimension(nfield,*)  yn,
integer, dimension(*)  ipkon,
integer, dimension(*)  inum,
integer, dimension(*)  kon,
character*8, dimension(*)  lakon,
integer  nfield,
integer  nk,
integer  ne,
integer, dimension(nfield,*)  nactdof 
)
21 !
22 ! interpolates basic degree of freedom nodal values
23 ! (displacements, temperatures) to 1d/2d nodal locations
24 !
25  implicit none
26 !
27  logical quadratic
28 !
29  character*8 lakon(*),lakonl
30 !
31  integer ipkon(*),inum(*),kon(*),ne,indexe,nfield,nk,i,j,k,l,
32  & node3(8,3),node6(3,6),node8(3,8),node2d,node3d,indexe2d,ne1d2d,
33  & node3m(8,3),iflag,nactdof(nfield,*),jmax
34 !
35  real*8 yn(nfield,*),ratioe(3)
36 !
37  include "gauss.f"
38 !
39  data node3 /1,4,8,5,12,20,16,17,9,11,15,13,
40  & 0,0,0,0,2,3,7,6,10,19,14,18/
41  data node3m /1,5,8,4,17,16,20,12,
42  & 0,0,0,0,0,0,0,0,
43  & 3,7,6,2,19,14,18,10/
44  data node6 /1,13,4,2,14,5,3,15,6,7,0,10,8,0,11,9,0,12/
45  data node8 /1,17,5,2,18,6,3,19,7,4,20,8,9,0,13,10,0,14,
46  & 11,0,15,12,0,16/
47  data ratioe /0.16666666666667d0,0.66666666666666d0,
48  & 0.16666666666667d0/
49  data iflag /2/
50 !
51 ! removing any results in 1d/2d nodes
52 !
53  ne1d2d=0
54 !
55  do i=1,ne
56 !
57  if(ipkon(i).lt.0) cycle
58  lakonl=lakon(i)
59  if((lakonl(7:7).eq.' ').or.(lakonl(7:7).eq.'I').or.
60  & (lakonl(1:1).ne.'C')) cycle
61  ne1d2d=1
62  indexe=ipkon(i)
63 c!
64 c! inactivating the 3d expansion nodes of 1d/2d elements
65 c!
66 c do j=1,20
67 c inum(kon(indexe+j))=0
68 c enddo
69 !
70  if((lakonl(4:5).eq.'15').or.(lakonl(4:4).eq.'6')) then
71  if(lakonl(4:5).eq.'15') then
72  indexe2d=indexe+15
73  jmax=6
74  else
75  indexe2d=indexe+6
76  jmax=3
77  endif
78  do j=1,jmax
79  node2d=kon(indexe2d+j)
80  inum(node2d)=0
81  do k=1,nfield
82  if(nactdof(k,node2d).le.0) yn(k,node2d)=0.d0
83  enddo
84  enddo
85  elseif(lakonl(7:7).eq.'B') then
86  if(lakonl(4:5).eq.'8I') then
87  indexe2d=indexe+11
88  jmax=2
89  elseif(lakonl(4:5).eq.'8R') then
90  indexe2d=indexe+8
91  jmax=2
92  elseif(lakonl(4:5).eq.'20') then
93  indexe2d=indexe+20
94  jmax=3
95  endif
96  do j=1,jmax
97  node2d=kon(indexe2d+j)
98  inum(node2d)=0
99  do k=1,nfield
100  if(nactdof(k,node2d).le.0) yn(k,node2d)=0.d0
101  enddo
102  enddo
103  else
104  if(lakonl(4:5).eq.'8I') then
105  indexe2d=indexe+11
106  jmax=4
107  elseif((lakonl(4:5).eq.'8R').or.(lakonl(4:5).eq.'8 ')) then
108  indexe2d=indexe+8
109  jmax=4
110  elseif(lakonl(4:5).eq.'20') then
111  indexe2d=indexe+20
112  jmax=8
113  endif
114  do j=1,jmax
115  node2d=kon(indexe2d+j)
116  inum(node2d)=0
117  do k=1,nfield
118  if(nactdof(k,node2d).le.0) yn(k,node2d)=0.d0
119  enddo
120  enddo
121  endif
122 !
123 ! inactivating the 3d expansion nodes of 1d/2d elements
124 !
125  do j=1,indexe2d-indexe
126  inum(kon(indexe+j))=0
127  enddo
128 !
129  enddo
130 !
131 ! if no 1d/2d elements return
132 !
133  if(ne1d2d.eq.0) return
134 !
135 ! interpolation of 3d results on 1d/2d nodes
136 !
137  do i=1,ne
138 !
139  if(ipkon(i).lt.0) cycle
140  lakonl=lakon(i)
141  if((lakonl(7:7).eq.' ').or.(lakonl(7:7).eq.'I').or.
142  & (lakonl(1:1).ne.'C')) cycle
143  indexe=ipkon(i)
144 !
145 ! check whether linear or quadratic element
146 !
147  if((lakonl(4:4).eq.'6').or.(lakonl(4:4).eq.'8')) then
148  quadratic=.false.
149  else
150  quadratic=.true.
151  endif
152 !
153  if((lakonl(4:5).eq.'15').or.(lakonl(4:4).eq.'6')) then
154  if(lakonl(4:5).eq.'15') then
155  indexe2d=indexe+15
156  jmax=6
157  else
158  indexe2d=indexe+6
159  jmax=3
160  endif
161  do j=1,jmax
162  node2d=kon(indexe2d+j)
163  inum(node2d)=inum(node2d)-1
164 !
165 ! taking the mean across the thickness
166 !
167  if((j.le.3).and.(quadratic)) then
168 !
169 ! end nodes: weights 1/6,2/3 and 1/6
170 !
171  do l=1,3
172  node3d=kon(indexe+node6(l,j))
173  do k=1,nfield
174  if(nactdof(k,node2d).le.0) yn(k,node2d)=
175  & yn(k,node2d)+yn(k,node3d)*ratioe(l)
176  enddo
177  enddo
178  else
179 !
180 ! middle nodes: weights 1/2,1/2
181 !
182  do l=1,3,2
183  node3d=kon(indexe+node6(l,j))
184  do k=1,nfield
185  if(nactdof(k,node2d).le.0) yn(k,node2d)=
186  & yn(k,node2d)+yn(k,node3d)/2.d0
187  enddo
188  enddo
189  endif
190  enddo
191  elseif(lakonl(7:7).eq.'B') then
192  if(lakonl(4:5).eq.'8I') then
193  indexe2d=indexe+11
194  jmax=2
195  elseif(lakonl(4:5).eq.'8R') then
196  indexe2d=indexe+8
197  jmax=2
198  elseif(lakonl(4:5).eq.'20') then
199  indexe2d=indexe+20
200  jmax=3
201  endif
202 !
203 ! mean values for beam elements
204 !
205  do j=1,jmax
206  node2d=kon(indexe2d+j)
207 !
208 ! mean value of vertex values
209 !
210  do l=1,4
211  inum(node2d)=inum(node2d)-1
212  if(quadratic) then
213  node3d=kon(indexe+node3(l,j))
214  else
215  node3d=kon(indexe+node3(l,2*j-1))
216  endif
217  do k=1,nfield
218  if(nactdof(k,node2d).le.0) yn(k,node2d)=
219  & yn(k,node2d)+yn(k,node3d)
220  enddo
221  enddo
222  enddo
223  else
224  if(lakonl(4:5).eq.'8I') then
225  indexe2d=indexe+11
226  jmax=4
227  elseif((lakonl(4:5).eq.'8R').or.(lakonl(4:5).eq.'8 ')) then
228  indexe2d=indexe+8
229  jmax=4
230  elseif(lakonl(4:5).eq.'20') then
231  indexe2d=indexe+20
232  jmax=8
233  endif
234  do j=1,jmax
235  node2d=kon(indexe2d+j)
236  inum(node2d)=inum(node2d)-1
237 !
238 ! taking the mean across the thickness
239 !
240  if((j.le.4).and.(quadratic)) then
241 !
242 ! end nodes: weights 1/6,2/3 and 1/6
243 !
244  do l=1,3
245  node3d=kon(indexe+node8(l,j))
246  do k=1,nfield
247  if(nactdof(k,node2d).le.0) yn(k,node2d)=
248  & yn(k,node2d)+yn(k,node3d)*ratioe(l)
249  enddo
250  enddo
251  else
252 !
253 ! middle nodes: weights 1/2,1/2
254 !
255  do l=1,3,2
256  node3d=kon(indexe+node8(l,j))
257  do k=1,nfield
258  if(nactdof(k,node2d).le.0) yn(k,node2d)=
259  & yn(k,node2d)+yn(k,node3d)/2.d0
260  enddo
261  enddo
262  endif
263  enddo
264  endif
265 !
266  enddo
267 !
268 ! taking the mean of nodal contributions coming from different
269 ! elements having the node in common
270 !
271  do i=1,nk
272  if(inum(i).lt.0) then
273  inum(i)=-inum(i)
274  do j=1,nfield
275  if(nactdof(j,i).le.0) yn(j,i)=yn(j,i)/inum(i)
276  enddo
277  endif
278  enddo
279 !
280  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)