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

Go to the source code of this file.

Functions/Subroutines

subroutine gen3dsurf (iponoel, inoel, iponoelmax, kon, ipkon, lakon, ne, iponor, knor, ipompc, nodempc, coefmpc, nmpc, nmpc_, mpcfree, ikmpc, ilmpc, labmpc, rig, ntrans, inotr, trab, nam, nk, nk_, co, nmethod, iperturb, nset, set, istartset, iendset, ialset, ikboun, ilboun, nboun, nboun_, nodeboun, ndirboun, xboun, iamboun, typeboun, mi, vold)
 

Function/Subroutine Documentation

◆ gen3dsurf()

subroutine gen3dsurf ( 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,
integer, dimension(*)  knor,
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(*)  rig,
integer  ntrans,
integer, dimension(2,*)  inotr,
real*8, dimension(7,*)  trab,
integer  nam,
integer  nk,
integer  nk_,
real*8, dimension(3,*)  co,
integer  nmethod,
integer  iperturb,
integer  nset,
character*81, dimension(*)  set,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
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(*)  mi,
real*8, dimension(0:mi(2),*)  vold 
)
25 !
26 ! connects nodes of 1-D and 2-D elements, belonging to nodal
27 ! surfaces, to the nodes of their expanded counterparts
28 !
29  implicit none
30 !
31  logical fixed
32 !
33  character*1 type,typeboun(*)
34  character*8 lakon(*)
35  character*20 labmpc(*),label
36  character*81 set(*)
37 !
38  integer iponoel(*),inoel(3,*),iponoelmax,kon(*),ipkon(*),ne,
39  & iponor(2,*),knor(*),ipompc(*),nodempc(3,*),nmpc,nmpc_,mpcfree,
40  & ikmpc(*),ilmpc(*),rig(*),ntrans,inotr(2,*),i,node,
41  & indexx,ielem,j,indexe,indexk,idir,nk,nk_,iamplitude,
42  & newnode,idof,id,mpcfreenew,k,nam,nmethod,iperturb,istartset(*),
43  & iendset(*),ialset(*),nset,ipos,l,idummy,ikboun(*),ilboun(*),
44  & nboun,nboun_,nodeboun(*),ndirboun(*),iamboun(*),mi(*)
45 !
46  real*8 coefmpc(*),trab(7,*),co(3,*),val,xboun(*),vold(0:mi(2),*)
47 !
48  label=' '
49  fixed=.false.
50 !
51  do i=1,nset
52  ipos=index(set(i),' ')
53 !
54 ! only nodal surfaces are involved
55 !
56  if(set(i)(ipos-1:ipos-1).ne.'S') cycle
57  do l=istartset(i),iendset(i)
58  node=ialset(l)
59  if(node.gt.iponoelmax) cycle
60  indexx=iponoel(node)
61  if(indexx.eq.0) cycle
62  ielem=inoel(1,indexx)
63  j=inoel(2,indexx)
64  indexe=ipkon(ielem)
65  indexk=iponor(2,indexe+j)
66 !
67  if(rig(node).eq.0) then
68 !
69 ! 2d element shell element: generate MPC's
70 !
71  if(lakon(ielem)(7:7).eq.'L') then
72  newnode=knor(indexk+1)
73  do idir=1,3
74  idof=8*(newnode-1)+idir
75  call nident(ikmpc,idof,nmpc,id)
76  if((id.le.0).or.(ikmpc(id).ne.idof)) then
77  nmpc=nmpc+1
78  if(nmpc.gt.nmpc_) then
79  write(*,*)
80  & '*ERROR in gen3dboun: increase nmpc_'
81  call exit(201)
82  endif
83  labmpc(nmpc)=' '
84  ipompc(nmpc)=mpcfree
85  do j=nmpc,id+2,-1
86  ikmpc(j)=ikmpc(j-1)
87  ilmpc(j)=ilmpc(j-1)
88  enddo
89  ikmpc(id+1)=idof
90  ilmpc(id+1)=nmpc
91  nodempc(1,mpcfree)=newnode
92  nodempc(2,mpcfree)=idir
93  coefmpc(mpcfree)=1.d0
94  mpcfree=nodempc(3,mpcfree)
95  if(mpcfree.eq.0) then
96  write(*,*)
97  & '*ERROR in gen3dboun: increase memmpc_'
98  call exit(201)
99  endif
100  nodempc(1,mpcfree)=knor(indexk+3)
101  nodempc(2,mpcfree)=idir
102  coefmpc(mpcfree)=1.d0
103  mpcfree=nodempc(3,mpcfree)
104  if(mpcfree.eq.0) then
105  write(*,*)
106  & '*ERROR in gen3dboun: increase memmpc_'
107  call exit(201)
108  endif
109  nodempc(1,mpcfree)=node
110  nodempc(2,mpcfree)=idir
111  coefmpc(mpcfree)=-2.d0
112  mpcfreenew=nodempc(3,mpcfree)
113  if(mpcfreenew.eq.0) then
114  write(*,*)
115  & '*ERROR in gen3dboun: increase memmpc_'
116  call exit(201)
117  endif
118  nodempc(3,mpcfree)=0
119  mpcfree=mpcfreenew
120  endif
121  enddo
122  elseif(lakon(ielem)(7:7).eq.'B') then
123 !
124 ! 1d beam element: generate MPC's
125 !
126  newnode=knor(indexk+1)
127  do idir=1,3
128  idof=8*(newnode-1)+idir
129  call nident(ikmpc,idof,nmpc,id)
130  if((id.le.0).or.(ikmpc(id).ne.idof)) then
131  nmpc=nmpc+1
132  if(nmpc.gt.nmpc_) then
133  write(*,*)
134  & '*ERROR in gen3dboun: increase nmpc_'
135  call exit(201)
136  endif
137  labmpc(nmpc)=' '
138  ipompc(nmpc)=mpcfree
139  do j=nmpc,id+2,-1
140  ikmpc(j)=ikmpc(j-1)
141  ilmpc(j)=ilmpc(j-1)
142  enddo
143  ikmpc(id+1)=idof
144  ilmpc(id+1)=nmpc
145  nodempc(1,mpcfree)=newnode
146  nodempc(2,mpcfree)=idir
147  coefmpc(mpcfree)=1.d0
148  mpcfree=nodempc(3,mpcfree)
149  if(mpcfree.eq.0) then
150  write(*,*)
151  & '*ERROR in gen3dboun: increase memmpc_'
152  call exit(201)
153  endif
154  do k=2,4
155  nodempc(1,mpcfree)=knor(indexk+k)
156  nodempc(2,mpcfree)=idir
157  coefmpc(mpcfree)=1.d0
158  mpcfree=nodempc(3,mpcfree)
159  if(mpcfree.eq.0) then
160  write(*,*)
161  & '*ERROR in gen3dboun: increase memmpc_'
162  call exit(201)
163  endif
164  enddo
165  nodempc(1,mpcfree)=node
166  nodempc(2,mpcfree)=idir
167  coefmpc(mpcfree)=-4.d0
168  mpcfreenew=nodempc(3,mpcfree)
169  if(mpcfreenew.eq.0) then
170  write(*,*)
171  & '*ERROR in gen3dboun: increase memmpc_'
172  call exit(201)
173  endif
174  nodempc(3,mpcfree)=0
175  mpcfree=mpcfreenew
176  endif
177 !
178  enddo
179  else
180 !
181 ! 2d plane stress, plane strain or axisymmetric
182 ! element: dependent node is replaced by new node in the middle
183 !
184 ! keeping the old node and generating an additional MPC leads
185 ! to problems since the old node is not restricted in the
186 ! z-direction, only the new node in the middle is. If the old
187 ! node is used subsequently in a contact spring element all
188 ! its DOFs are used and the unrestricted z-DOF leads to a
189 ! singular equation system
190 !
191 c write(*,*) ialset(l),' replaced by ',knor(indexk+2)
192 c co(1,knor(indexk+2))=co(1,ialset(l))
193 c co(2,knor(indexk+2))=co(2,ialset(l))
194 c co(3,knor(indexk+2))=co(3,ialset(l))
195 c ialset(l)=knor(indexk+2)
196 !
197 ! 2d plane stress, plane strain or axisymmetric
198 ! element: MPC in all but z-direction
199 !
200  newnode=knor(indexk+2)
201  do idir=1,2
202  idof=8*(newnode-1)+idir
203  call nident(ikmpc,idof,nmpc,id)
204  if(((id.le.0).or.(ikmpc(id).ne.idof)).and.
205  & (idir.ne.3)) then
206  nmpc=nmpc+1
207  if(nmpc.gt.nmpc_) then
208  write(*,*)
209  & '*ERROR in gen3dmpc: increase nmpc_'
210  call exit(201)
211  endif
212  labmpc(nmpc)=' '
213  ipompc(nmpc)=mpcfree
214  do j=nmpc,id+2,-1
215  ikmpc(j)=ikmpc(j-1)
216  ilmpc(j)=ilmpc(j-1)
217  enddo
218  ikmpc(id+1)=idof
219  ilmpc(id+1)=nmpc
220  nodempc(1,mpcfree)=newnode
221  nodempc(2,mpcfree)=idir
222  coefmpc(mpcfree)=1.d0
223  mpcfree=nodempc(3,mpcfree)
224  if(mpcfree.eq.0) then
225  write(*,*)
226  & '*ERROR in gen3dmpc: increase memmpc_'
227  call exit(201)
228  endif
229  nodempc(1,mpcfree)=node
230  nodempc(2,mpcfree)=idir
231  coefmpc(mpcfree)=-1.d0
232  mpcfreenew=nodempc(3,mpcfree)
233  if(mpcfreenew.eq.0) then
234  write(*,*)
235  & '*ERROR in gen3dmpc: increase memmpc_'
236  call exit(201)
237  endif
238  nodempc(3,mpcfree)=0
239  mpcfree=mpcfreenew
240  endif
241  enddo
242 !
243 ! fixing the original node in z-direction
244 !
245  val=0.d0
246  k=3
247  if(nam.gt.0) iamplitude=0
248  type='M'
249  call bounadd(node,k,k,val,nodeboun,
250  & ndirboun,xboun,nboun,nboun_,iamboun,
251  & iamplitude,nam,ipompc,nodempc,coefmpc,
252  & nmpc,nmpc_,mpcfree,inotr,trab,ntrans,
253  & ikboun,ilboun,ikmpc,ilmpc,co,nk,nk_,
254  & labmpc,type,typeboun,nmethod,iperturb,
255  & fixed,vold,idummy,mi,label)
256 !
257  endif
258  endif
259  enddo
260  enddo
261 !
262  return
subroutine nident(x, px, n, id)
Definition: nident.f:26
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
Hosted by OpenAircraft.com, (Michigan UAV, LLC)