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

Go to the source code of this file.

Functions/Subroutines

subroutine straightmpc (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

◆ straightmpc()

subroutine straightmpc ( 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
38 !
39  real*8 coefmpc(3,*),co(3,*),dd,dmax,xboun(*)
40 !
41  save nodea,nodeb,jmax
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 k=1,3
50  dd=abs((co(k,nodea)-co(k,nodeb)))
51  if(dd.gt.dmax) then
52  dmax=dd
53  jmax=k
54  endif
55  enddo
56  return
57  endif
58 !
59  nk=nk+1
60  if(nk.gt.nk_) then
61  write(*,*) '*ERROR in straightmpc: increase nk_'
62  call exit(201)
63  endif
64  do j=1,3
65  if(j.eq.jmax) cycle
66  idof=8*(node-1)+j
67  call nident(ikmpc,idof,nmpc,id)
68  if(id.gt.0) then
69  if(ikmpc(id).eq.idof) then
70  write(*,*) '*WARNING in straightmpc: DOF for node ',node
71  write(*,*) ' in direction ',j,' has been used'
72  write(*,*) ' on the dependent side of another MPC'
73  write(*,*) ' STRAIGHT constraint cannot be applied'
74  cycle
75  endif
76  endif
77  nmpc=nmpc+1
78  if(nmpc.gt.nmpc_) then
79  write(*,*) '*ERROR in straightmpc: increase nmpc_'
80  call exit(201)
81  endif
82 !
83  ipompc(nmpc)=mpcfree
84  labmpc(nmpc)='STRAIGHT '
85 !
86  do l=nmpc,id+2,-1
87  ikmpc(l)=ikmpc(l-1)
88  ilmpc(l)=ilmpc(l-1)
89  enddo
90  ikmpc(id+1)=idof
91  ilmpc(id+1)=nmpc
92 !
93  nodempc(1,mpcfree)=node
94  nodempc(2,mpcfree)=j
95  mpcfree=nodempc(3,mpcfree)
96  nodempc(1,mpcfree)=node
97  nodempc(2,mpcfree)=jmax
98  mpcfree=nodempc(3,mpcfree)
99  nodempc(1,mpcfree)=nodea
100  nodempc(2,mpcfree)=j
101  mpcfree=nodempc(3,mpcfree)
102  nodempc(1,mpcfree)=nodea
103  nodempc(2,mpcfree)=jmax
104  mpcfree=nodempc(3,mpcfree)
105  nodempc(1,mpcfree)=nodeb
106  nodempc(2,mpcfree)=j
107  mpcfree=nodempc(3,mpcfree)
108  nodempc(1,mpcfree)=nodeb
109  nodempc(2,mpcfree)=jmax
110  mpcfree=nodempc(3,mpcfree)
111  nodempc(1,mpcfree)=nk
112  nodempc(2,mpcfree)=j
113  mpcfreeold=mpcfree
114  mpcfree=nodempc(3,mpcfree)
115  nodempc(3,mpcfreeold)=0
116  idof=8*(nk-1)+j
117  call nident(ikboun,idof,nboun,id)
118  nboun=nboun+1
119  if(nboun.gt.nboun_) then
120  write(*,*) '*ERROR in straightmpc: increase nboun_'
121  call exit(201)
122  endif
123  nodeboun(nboun)=nk
124  ndirboun(nboun)=j
125  typeboun(nboun)='U'
126  do l=nboun,id+2,-1
127  ikboun(l)=ikboun(l-1)
128  ilboun(l)=ilboun(l-1)
129  enddo
130  ikboun(id+1)=idof
131  ilboun(id+1)=nboun
132  enddo
133 !
134  return
subroutine nident(x, px, n, id)
Definition: nident.f:26
Hosted by OpenAircraft.com, (Michigan UAV, LLC)