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

Go to the source code of this file.

Functions/Subroutines

subroutine beammpc (ipompc, nodempc, labmpc, nmpc, nmpc_, mpcfree, ikmpc, ilmpc, nk, nk_, nodeboun, ndirboun, ikboun, ilboun, nboun, nboun_, inode, node, co, typeboun)
 

Function/Subroutine Documentation

◆ beammpc()

subroutine beammpc ( integer, dimension(*)  ipompc,
integer, dimension(3,*)  nodempc,
character*20, dimension(*)  labmpc,
integer  nmpc,
integer  nmpc_,
integer  mpcfree,
integer, dimension(*)  ikmpc,
integer, dimension(*)  ilmpc,
integer  nk,
integer  nk_,
integer, dimension(*)  nodeboun,
integer, dimension(*)  ndirboun,
integer, dimension(*)  ikboun,
integer, dimension(*)  ilboun,
integer  nboun,
integer  nboun_,
integer  inode,
integer  node,
real*8, dimension(3,*)  co,
character*1, dimension(*)  typeboun 
)
22 !
23 ! generates an MPC for two nodes staying on a constant distance
24 ! from each other. The parameter inode indicates how many
25 ! times the present routine was called within the same *MPC
26 ! definition. For inode=1 "node" is node a, for inode=2 "node"
27 ! is node b.
28 !
29  implicit none
30 !
31  character*1 typeboun(*)
32  character*20 labmpc(*)
33 !
34  integer ipompc(*),nodempc(3,*),nmpc,nmpc_,mpcfree,nk,nk_,ikmpc(*),
35  & ilmpc(*),node,id,mpcfreeold,j,idof,l,nodeboun(*),nodea,nodeb,
36  & ndirboun(*),ikboun(*),ilboun(*),nboun,nboun_,inode,jmax,k,
37  & m
38 !
39  real*8 co(3,*),dd,dmax
40 !
41  save nodea
42 !
43  if(inode.eq.1) then
44  nodea=node
45  return
46  elseif(inode.eq.2) then
47  nodeb=node
48  dmax=0.d0
49  do j=1,3
50  dd=dabs(co(j,nodea)-co(j,nodeb))
51  if(dd.gt.dmax) then
52  dmax=dd
53  jmax=j
54  endif
55  enddo
56  endif
57 !
58  nk=nk+1
59  if(nk.gt.nk_) then
60  write(*,*) '*ERROR in beammpc: increase nk_'
61  call exit(201)
62  endif
63 !
64  j=jmax
65  k=j+1
66  if(k.gt.3) k=1
67  l=k+1
68  if(l.gt.3) l=1
69 !
70  idof=8*(nodea-1)+j
71  call nident(ikmpc,idof,nmpc,id)
72  if(id.gt.0) then
73  if(ikmpc(id).eq.idof) then
74  write(*,*) '*WARNING in beammpc: DOF for node ',node
75  write(*,*) ' in direction ',j,' has been used'
76  write(*,*) ' on the dependent side of another MPC'
77  write(*,*) ' PLANE constraint cannot be applied'
78  return
79  endif
80  endif
81  nmpc=nmpc+1
82  if(nmpc.gt.nmpc_) then
83  write(*,*) '*ERROR in beammpc: increase nmpc_'
84  call exit(201)
85  endif
86 !
87  ipompc(nmpc)=mpcfree
88  labmpc(nmpc)='BEAM '
89 !
90  do m=nmpc,id+2,-1
91  ikmpc(m)=ikmpc(m-1)
92  ilmpc(m)=ilmpc(m-1)
93  enddo
94  ikmpc(id+1)=idof
95  ilmpc(id+1)=nmpc
96 !
97  nodempc(1,mpcfree)=nodea
98  nodempc(2,mpcfree)=j
99  mpcfree=nodempc(3,mpcfree)
100  nodempc(1,mpcfree)=nodea
101  nodempc(2,mpcfree)=k
102  mpcfree=nodempc(3,mpcfree)
103  nodempc(1,mpcfree)=nodea
104  nodempc(2,mpcfree)=l
105  mpcfree=nodempc(3,mpcfree)
106  nodempc(1,mpcfree)=nodeb
107  nodempc(2,mpcfree)=j
108  mpcfree=nodempc(3,mpcfree)
109  nodempc(1,mpcfree)=nodeb
110  nodempc(2,mpcfree)=k
111  mpcfree=nodempc(3,mpcfree)
112  nodempc(1,mpcfree)=nodeb
113  nodempc(2,mpcfree)=l
114  mpcfree=nodempc(3,mpcfree)
115  nodempc(1,mpcfree)=nk
116  nodempc(2,mpcfree)=j
117  mpcfreeold=mpcfree
118  mpcfree=nodempc(3,mpcfree)
119  nodempc(3,mpcfreeold)=0
120  idof=8*(nk-1)+j
121  call nident(ikboun,idof,nboun,id)
122  nboun=nboun+1
123  if(nboun.gt.nboun_) then
124  write(*,*) '*ERROR in beammpc: increase nboun_'
125  call exit(201)
126  endif
127  nodeboun(nboun)=nk
128  ndirboun(nboun)=j
129  typeboun(nboun)='U'
130  do m=nboun,id+2,-1
131  ikboun(m)=ikboun(m-1)
132  ilboun(m)=ilboun(m-1)
133  enddo
134  ikboun(id+1)=idof
135  ilboun(id+1)=nboun
136 !
137  return
subroutine nident(x, px, n, id)
Definition: nident.f:26
Hosted by OpenAircraft.com, (Michigan UAV, LLC)