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

Go to the source code of this file.

Functions/Subroutines

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

Function/Subroutine Documentation

◆ planempc()

subroutine planempc ( integer, dimension(*)  ipompc,
integer, dimension(3,*)  nodempc,
real*8, dimension(3,*)  coefmpc,
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_,
real*8, dimension(*)  xboun,
integer  inode,
integer  node,
real*8, dimension(3,*)  co,
character*1, dimension(*)  typeboun 
)
22 !
23 ! generates MPC's for nodes staying on a straight line defined
24 ! by two nodes a and b. 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. Starting with inode=3 MPC's are defined.
28 !
29  implicit none
30 !
31  character*1 typeboun(*)
32  character*20 labmpc(*)
33 !
34  integer ipompc(*),nodempc(3,*),nmpc,nmpc_,mpcfree,nk,nk_,
35  & ikmpc(*),
36  & ilmpc(*),node,id,mpcfreeold,j,idof,l,nodeboun(*),nodea,nodeb,
37  & ndirboun(*),ikboun(*),ilboun(*),nboun,nboun_,inode,jmax,k,nodec,
38  & m
39 !
40  real*8 coefmpc(3,*),co(3,*),dd,dmax,pac(3),pbc(3),xboun(*)
41 !
42  save nodea,nodeb,nodec,jmax
43 !
44  if(inode.eq.1) then
45  nodea=node
46  return
47  elseif(inode.eq.2) then
48  nodeb=node
49  return
50  elseif(inode.eq.3) then
51  nodec=node
52  do j=1,3
53  pac(j)=co(j,nodea)-co(j,nodec)
54  pbc(j)=co(j,nodeb)-co(j,nodec)
55  enddo
56  dmax=abs(pac(2)*pbc(3)-pac(3)*pbc(2))
57  jmax=1
58  dd=abs(pac(1)*pbc(3)-pac(3)*pbc(1))
59  if(dd.gt.dmax) then
60  dmax=dd
61  jmax=2
62  endif
63  dd=abs(pac(1)*pbc(2)-pac(2)*pbc(1))
64  if(dd.gt.dmax) then
65  dmax=dd
66  jmax=3
67  endif
68  return
69  endif
70 !
71  nk=nk+1
72  if(nk.gt.nk_) then
73  write(*,*) '*ERROR in planempc: increase nk_'
74  call exit(201)
75  endif
76 !
77  j=jmax
78  k=j+1
79  if(k.gt.3) k=1
80  l=k+1
81  if(l.gt.3) l=1
82 !
83  idof=8*(node-1)+j
84  call nident(ikmpc,idof,nmpc,id)
85  if(id.gt.0) then
86  if(ikmpc(id).eq.idof) then
87  write(*,*) '*WARNING in planempc: DOF for node ',node
88  write(*,*) ' in direction ',j,' has been used'
89  write(*,*) ' on the dependent side of another MPC'
90  write(*,*) ' PLANE constraint cannot be applied'
91  return
92  endif
93  endif
94  nmpc=nmpc+1
95  if(nmpc.gt.nmpc_) then
96  write(*,*) '*ERROR in planempc: increase nmpc_'
97  call exit(201)
98  endif
99 !
100  ipompc(nmpc)=mpcfree
101  labmpc(nmpc)='PLANE '
102 !
103  do m=nmpc,id+2,-1
104  ikmpc(m)=ikmpc(m-1)
105  ilmpc(m)=ilmpc(m-1)
106  enddo
107  ikmpc(id+1)=idof
108  ilmpc(id+1)=nmpc
109 !
110  nodempc(1,mpcfree)=node
111  nodempc(2,mpcfree)=j
112  mpcfree=nodempc(3,mpcfree)
113  nodempc(1,mpcfree)=node
114  nodempc(2,mpcfree)=k
115  mpcfree=nodempc(3,mpcfree)
116  nodempc(1,mpcfree)=node
117  nodempc(2,mpcfree)=l
118  mpcfree=nodempc(3,mpcfree)
119  nodempc(1,mpcfree)=nodea
120  nodempc(2,mpcfree)=j
121  mpcfree=nodempc(3,mpcfree)
122  nodempc(1,mpcfree)=nodea
123  nodempc(2,mpcfree)=k
124  mpcfree=nodempc(3,mpcfree)
125  nodempc(1,mpcfree)=nodea
126  nodempc(2,mpcfree)=l
127  mpcfree=nodempc(3,mpcfree)
128  nodempc(1,mpcfree)=nodeb
129  nodempc(2,mpcfree)=j
130  mpcfree=nodempc(3,mpcfree)
131  nodempc(1,mpcfree)=nodeb
132  nodempc(2,mpcfree)=k
133  mpcfree=nodempc(3,mpcfree)
134  nodempc(1,mpcfree)=nodeb
135  nodempc(2,mpcfree)=l
136  mpcfree=nodempc(3,mpcfree)
137  nodempc(1,mpcfree)=nodec
138  nodempc(2,mpcfree)=j
139  mpcfree=nodempc(3,mpcfree)
140  nodempc(1,mpcfree)=nodec
141  nodempc(2,mpcfree)=k
142  mpcfree=nodempc(3,mpcfree)
143  nodempc(1,mpcfree)=nodec
144  nodempc(2,mpcfree)=l
145  mpcfree=nodempc(3,mpcfree)
146  nodempc(1,mpcfree)=nk
147  nodempc(2,mpcfree)=j
148  mpcfreeold=mpcfree
149  mpcfree=nodempc(3,mpcfree)
150  nodempc(3,mpcfreeold)=0
151  idof=8*(nk-1)+j
152  call nident(ikboun,idof,nboun,id)
153  nboun=nboun+1
154  if(nboun.gt.nboun_) then
155  write(*,*) '*ERROR in planempc: increase nboun_'
156  call exit(201)
157  endif
158  nodeboun(nboun)=nk
159  ndirboun(nboun)=j
160  typeboun(nboun)='U'
161  do m=nboun,id+2,-1
162  ikboun(m)=ikboun(m-1)
163  ilboun(m)=ilboun(m-1)
164  enddo
165  ikboun(id+1)=idof
166  ilboun(id+1)=nboun
167 !
168  return
subroutine nident(x, px, n, id)
Definition: nident.f:26
Hosted by OpenAircraft.com, (Michigan UAV, LLC)