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

Go to the source code of this file.

Functions/Subroutines

subroutine gen3dmpc (ipompc, nodempc, coefmpc, nmpc, nmpc_, mpcfree, ikmpc, ilmpc, labmpc, iponoel, inoel, iponoelmax, kon, ipkon, lakon, ne, iponor, xnor, knor, rig)
 

Function/Subroutine Documentation

◆ gen3dmpc()

subroutine gen3dmpc ( 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, dimension(*)  iponoel,
integer, dimension(3,*)  inoel,
integer  iponoelmax,
integer, dimension(*)  kon,
integer, dimension(*)  ipkon,
character*8, dimension(*)  lakon,
integer  ne,
integer, dimension(2,*)  iponor,
real*8, dimension(*)  xnor,
integer, dimension(*)  knor,
integer, dimension(*)  rig 
)
22 !
23 ! connects nodes of 1-D and 2-D elements, for which MPC's were
24 ! defined, to the nodes of their expanded counterparts
25 !
26  implicit none
27 !
28  logical dependent
29 !
30  character*8 lakon(*)
31  character*20 labmpc(*)
32 !
33  integer ipompc(*),nodempc(3,*),nmpc,nmpc_,mpcfree,ikmpc(*),
34  & ilmpc(*),iponoel(*),inoel(3,*),iponoelmax,kon(*),ipkon(*),
35  & ne,iponor(2,*),knor(*),rig(*),i,index1,node,index2,ielem,
36  & indexe,j,indexk,newnode,idir,idof,id,mpcfreenew,k,idofold,
37  & idofnew,idold,idnew
38 !
39  real*8 coefmpc(*),xnor(*)
40 !
41  do i=1,nmpc
42  index1=ipompc(i)
43  dependent=.true.
44  do
45  node=nodempc(1,index1)
46  if(node.le.iponoelmax) then
47  if(rig(node).ne.0) then
48  if(nodempc(2,index1).gt.3) then
49  if(rig(node).lt.0) then
50  write(*,*) '*ERROR in gen3dmpc: in node ',node
51  write(*,*) ' a rotational DOF is constrained'
52  write(*,*) ' by a SPC; however, the elements'
53  write(*,*) ' to which this node belongs do not'
54  write(*,*) ' have rotational DOFs'
55  call exit(201)
56  endif
57  nodempc(1,index1)=rig(node)
58  nodempc(2,index1)=nodempc(2,index1)-3
59 !
60 ! adapting ikmpc and ilmpc (only for the dependent
61 ! term)
62 !
63  if(dependent) then
64  idofold=8*(node-1)+nodempc(2,index1)+3
65  call nident(ikmpc,idofold,nmpc,idold)
66  idofnew=8*(rig(node)-1)+nodempc(2,index1)
67  call nident(ikmpc,idofnew,nmpc,idnew)
68  if(idold.le.idnew) then
69  do j=idold,idnew-1
70  ikmpc(j)=ikmpc(j+1)
71  ilmpc(j)=ilmpc(j+1)
72  enddo
73  else
74  do j=idold,idnew+1,-1
75  ikmpc(j)=ikmpc(j-1)
76  ilmpc(j)=ilmpc(j-1)
77  enddo
78  endif
79  ikmpc(idnew)=idofnew
80  ilmpc(idnew)=i
81  endif
82 !
83  endif
84  else
85  index2=iponoel(node)
86 !
87 ! check for nodes not belonging to 1d or 2d elements
88 !
89  if(index2.eq.0) then
90  index1=nodempc(3,index1)
91  dependent=.false.
92  if(index1.eq.0) exit
93  cycle
94  endif
95 c
96  ielem=inoel(1,index2)
97  indexe=ipkon(ielem)
98  j=inoel(2,index2)
99  indexk=iponor(2,indexe+j)
100 !
101 ! 2d element shell element
102 !
103  if(lakon(ielem)(7:7).eq.'L') then
104  newnode=knor(indexk+1)
105  idir=nodempc(2,index1)
106  idof=8*(newnode-1)+idir
107  call nident(ikmpc,idof,nmpc,id)
108  if((id.le.0).or.(ikmpc(id).ne.idof)) then
109  nmpc=nmpc+1
110  if(nmpc.gt.nmpc_) then
111  write(*,*)
112  & '*ERROR in gen3dmpc: increase nmpc_'
113  call exit(201)
114  endif
115  labmpc(nmpc)=' '
116  ipompc(nmpc)=mpcfree
117  do j=nmpc,id+2,-1
118  ikmpc(j)=ikmpc(j-1)
119  ilmpc(j)=ilmpc(j-1)
120  enddo
121  ikmpc(id+1)=idof
122  ilmpc(id+1)=nmpc
123  nodempc(1,mpcfree)=newnode
124  nodempc(2,mpcfree)=idir
125  coefmpc(mpcfree)=1.d0
126  mpcfree=nodempc(3,mpcfree)
127  if(mpcfree.eq.0) then
128  write(*,*)
129  & '*ERROR in gen3dmpc: increase memmpc_'
130  call exit(201)
131  endif
132  nodempc(1,mpcfree)=knor(indexk+3)
133  nodempc(2,mpcfree)=idir
134  coefmpc(mpcfree)=1.d0
135  mpcfree=nodempc(3,mpcfree)
136  if(mpcfree.eq.0) then
137  write(*,*)
138  & '*ERROR in gen3dmpc: increase memmpc_'
139  call exit(201)
140  endif
141  nodempc(1,mpcfree)=node
142  nodempc(2,mpcfree)=idir
143  coefmpc(mpcfree)=-2.d0
144  mpcfreenew=nodempc(3,mpcfree)
145  if(mpcfreenew.eq.0) then
146  write(*,*)
147  & '*ERROR in gen3dmpc: increase memmpc_'
148  call exit(201)
149  endif
150  nodempc(3,mpcfree)=0
151  mpcfree=mpcfreenew
152  endif
153 !
154 ! fixing the temperature degrees of freedom
155 !
156  if(idir.eq.0) then
157 !
158 ! t(n_3)=t(n)
159 !
160  newnode=knor(indexk+3)
161  idof=8*(newnode-1)+idir
162  call nident(ikmpc,idof,nmpc,id)
163  if((id.le.0).or.(ikmpc(id).ne.idof)) then
164  nmpc=nmpc+1
165  if(nmpc.gt.nmpc_) then
166  write(*,*)
167  & '*ERROR in gen3dboun: increase nmpc_'
168  call exit(201)
169  endif
170  labmpc(nmpc)=' '
171  ipompc(nmpc)=mpcfree
172  do j=nmpc,id+2,-1
173  ikmpc(j)=ikmpc(j-1)
174  ilmpc(j)=ilmpc(j-1)
175  enddo
176  ikmpc(id+1)=idof
177  ilmpc(id+1)=nmpc
178  nodempc(1,mpcfree)=newnode
179  nodempc(2,mpcfree)=idir
180  coefmpc(mpcfree)=1.d0
181  mpcfree=nodempc(3,mpcfree)
182  if(mpcfree.eq.0) then
183  write(*,*)
184  & '*ERROR in gen3dboun: increase memmpc_'
185  call exit(201)
186  endif
187  nodempc(1,mpcfree)=node
188  nodempc(2,mpcfree)=idir
189  coefmpc(mpcfree)=-1.d0
190  mpcfreenew=nodempc(3,mpcfree)
191  if(mpcfreenew.eq.0) then
192  write(*,*)
193  & '*ERROR in gen3dboun: increase memmpc_'
194  call exit(201)
195  endif
196  nodempc(3,mpcfree)=0
197  mpcfree=mpcfreenew
198  endif
199  endif
200  elseif(lakon(ielem)(7:7).eq.'B') then
201 !
202 ! 1d beam element
203 !
204  newnode=knor(indexk+1)
205  idir=nodempc(2,index1)
206  idof=8*(newnode-1)+idir
207  call nident(ikmpc,idof,nmpc,id)
208  if((id.le.0).or.(ikmpc(id).ne.idof)) then
209  nmpc=nmpc+1
210  if(nmpc.gt.nmpc_) then
211  write(*,*)
212  & '*ERROR in gen3dmpc: increase nmpc_'
213  call exit(201)
214  endif
215  labmpc(nmpc)=' '
216  ipompc(nmpc)=mpcfree
217  do j=nmpc,id+2,-1
218  ikmpc(j)=ikmpc(j-1)
219  ilmpc(j)=ilmpc(j-1)
220  enddo
221  ikmpc(id+1)=idof
222  ilmpc(id+1)=nmpc
223  nodempc(1,mpcfree)=newnode
224  nodempc(2,mpcfree)=idir
225  coefmpc(mpcfree)=1.d0
226  mpcfree=nodempc(3,mpcfree)
227  if(mpcfree.eq.0) then
228  write(*,*)
229  & '*ERROR in gen3dmpc: increase memmpc_'
230  call exit(201)
231  endif
232  do k=2,4
233  nodempc(1,mpcfree)=knor(indexk+k)
234  nodempc(2,mpcfree)=idir
235  coefmpc(mpcfree)=1.d0
236  mpcfree=nodempc(3,mpcfree)
237  if(mpcfree.eq.0) then
238  write(*,*)
239  & '*ERROR in gen3dmpc: increase memmpc_'
240  call exit(201)
241  endif
242  enddo
243  nodempc(1,mpcfree)=node
244  nodempc(2,mpcfree)=idir
245  coefmpc(mpcfree)=-4.d0
246  mpcfreenew=nodempc(3,mpcfree)
247  if(mpcfreenew.eq.0) then
248  write(*,*)
249  & '*ERROR in gen3dmpc: increase memmpc_'
250  call exit(201)
251  endif
252  nodempc(3,mpcfree)=0
253  mpcfree=mpcfreenew
254  endif
255 !
256 ! fixing the temperature degrees of freedom
257 !
258  if(idir.eq.0) then
259  do k=2,4
260 !
261 ! t(n_k)=t(n), k=2,4
262 !
263  newnode=knor(indexk+k)
264  idof=8*(newnode-1)+idir
265  call nident(ikmpc,idof,nmpc,id)
266  if((id.le.0).or.(ikmpc(id).ne.idof)) then
267  nmpc=nmpc+1
268  if(nmpc.gt.nmpc_) then
269  write(*,*)
270  & '*ERROR in gen3dboun: increase nmpc_'
271  call exit(201)
272  endif
273  labmpc(nmpc)=' '
274  ipompc(nmpc)=mpcfree
275  do j=nmpc,id+2,-1
276  ikmpc(j)=ikmpc(j-1)
277  ilmpc(j)=ilmpc(j-1)
278  enddo
279  ikmpc(id+1)=idof
280  ilmpc(id+1)=nmpc
281  nodempc(1,mpcfree)=newnode
282  nodempc(2,mpcfree)=idir
283  coefmpc(mpcfree)=1.d0
284  mpcfree=nodempc(3,mpcfree)
285  if(mpcfree.eq.0) then
286  write(*,*)
287  & '*ERROR in gen3dboun: increase memmpc_'
288  call exit(201)
289  endif
290  nodempc(1,mpcfree)=node
291  nodempc(2,mpcfree)=idir
292  coefmpc(mpcfree)=-1.d0
293  mpcfreenew=nodempc(3,mpcfree)
294  if(mpcfreenew.eq.0) then
295  write(*,*)
296  & '*ERROR in gen3dboun: increase memmpc_'
297  call exit(201)
298  endif
299  nodempc(3,mpcfree)=0
300  mpcfree=mpcfreenew
301  endif
302  enddo
303  endif
304  else
305 !
306 ! 2d plane stress, plane strain or axisymmetric
307 ! element
308 !
309  newnode=knor(indexk+2)
310  idir=nodempc(2,index1)
311  idof=8*(newnode-1)+idir
312  call nident(ikmpc,idof,nmpc,id)
313  if(((id.le.0).or.(ikmpc(id).ne.idof)).and.
314  & (idir.ne.3)) then
315  nmpc=nmpc+1
316  if(nmpc.gt.nmpc_) then
317  write(*,*)
318  & '*ERROR in gen3dmpc: increase nmpc_'
319  call exit(201)
320  endif
321  labmpc(nmpc)=' '
322  ipompc(nmpc)=mpcfree
323  do j=nmpc,id+2,-1
324  ikmpc(j)=ikmpc(j-1)
325  ilmpc(j)=ilmpc(j-1)
326  enddo
327  ikmpc(id+1)=idof
328  ilmpc(id+1)=nmpc
329  nodempc(1,mpcfree)=newnode
330  nodempc(2,mpcfree)=idir
331  coefmpc(mpcfree)=1.d0
332  mpcfree=nodempc(3,mpcfree)
333  if(mpcfree.eq.0) then
334  write(*,*)
335  & '*ERROR in gen3dmpc: increase memmpc_'
336  call exit(201)
337  endif
338  nodempc(1,mpcfree)=node
339  nodempc(2,mpcfree)=idir
340  coefmpc(mpcfree)=-1.d0
341  mpcfreenew=nodempc(3,mpcfree)
342  if(mpcfreenew.eq.0) then
343  write(*,*)
344  & '*ERROR in gen3dmpc: increase memmpc_'
345  call exit(201)
346  endif
347  nodempc(3,mpcfree)=0
348  mpcfree=mpcfreenew
349  endif
350  endif
351  endif
352  endif
353  index1=nodempc(3,index1)
354  dependent=.false.
355  if(index1.eq.0) exit
356  enddo
357  enddo
358 !
359  return
subroutine nident(x, px, n, id)
Definition: nident.f:26
Hosted by OpenAircraft.com, (Michigan UAV, LLC)